1 // Copyright (C) 2007-2012 CEA/DEN, EDF R&D, OPEN CASCADE
3 // Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License.
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
22 #include "SMESH_ControlsDef.hxx"
24 #include "SMDS_BallElement.hxx"
25 #include "SMDS_Iterator.hxx"
26 #include "SMDS_Mesh.hxx"
27 #include "SMDS_MeshElement.hxx"
28 #include "SMDS_MeshNode.hxx"
29 #include "SMDS_QuadraticEdge.hxx"
30 #include "SMDS_QuadraticFaceOfNodes.hxx"
31 #include "SMDS_VolumeTool.hxx"
32 #include "SMESHDS_GroupBase.hxx"
33 #include "SMESHDS_Mesh.hxx"
34 #include "SMESH_OctreeNode.hxx"
36 #include <BRepAdaptor_Surface.hxx>
37 #include <BRepClass_FaceClassifier.hxx>
38 #include <BRep_Tool.hxx>
39 #include <Geom_CylindricalSurface.hxx>
40 #include <Geom_Plane.hxx>
41 #include <Geom_Surface.hxx>
42 #include <Precision.hxx>
43 #include <TColStd_MapIteratorOfMapOfInteger.hxx>
44 #include <TColStd_MapOfInteger.hxx>
45 #include <TColStd_SequenceOfAsciiString.hxx>
46 #include <TColgp_Array1OfXYZ.hxx>
49 #include <TopoDS_Edge.hxx>
50 #include <TopoDS_Face.hxx>
51 #include <TopoDS_Iterator.hxx>
52 #include <TopoDS_Shape.hxx>
53 #include <TopoDS_Vertex.hxx>
55 #include <gp_Cylinder.hxx>
62 #include <vtkMeshQuality.h>
74 const double theEps = 1e-100;
75 const double theInf = 1e+100;
77 inline gp_XYZ gpXYZ(const SMDS_MeshNode* aNode )
79 return gp_XYZ(aNode->X(), aNode->Y(), aNode->Z() );
82 inline double getAngle( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
84 gp_Vec v1( P1 - P2 ), v2( P3 - P2 );
86 return v1.Magnitude() < gp::Resolution() ||
87 v2.Magnitude() < gp::Resolution() ? 0 : v1.Angle( v2 );
90 inline double getArea( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
92 gp_Vec aVec1( P2 - P1 );
93 gp_Vec aVec2( P3 - P1 );
94 return ( aVec1 ^ aVec2 ).Magnitude() * 0.5;
97 inline double getArea( const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& P3 )
99 return getArea( P1.XYZ(), P2.XYZ(), P3.XYZ() );
104 inline double getDistance( const gp_XYZ& P1, const gp_XYZ& P2 )
106 double aDist = gp_Pnt( P1 ).Distance( gp_Pnt( P2 ) );
110 int getNbMultiConnection( const SMDS_Mesh* theMesh, const int theId )
115 const SMDS_MeshElement* anEdge = theMesh->FindElement( theId );
116 if ( anEdge == 0 || anEdge->GetType() != SMDSAbs_Edge/* || anEdge->NbNodes() != 2 */)
119 // for each pair of nodes in anEdge (there are 2 pairs in a quadratic edge)
120 // count elements containing both nodes of the pair.
121 // Note that there may be such cases for a quadratic edge (a horizontal line):
126 // +-----+------+ +-----+------+
129 // result sould be 2 in both cases
131 int aResult0 = 0, aResult1 = 0;
132 // last node, it is a medium one in a quadratic edge
133 const SMDS_MeshNode* aLastNode = anEdge->GetNode( anEdge->NbNodes() - 1 );
134 const SMDS_MeshNode* aNode0 = anEdge->GetNode( 0 );
135 const SMDS_MeshNode* aNode1 = anEdge->GetNode( 1 );
136 if ( aNode1 == aLastNode ) aNode1 = 0;
138 SMDS_ElemIteratorPtr anElemIter = aLastNode->GetInverseElementIterator();
139 while( anElemIter->more() ) {
140 const SMDS_MeshElement* anElem = anElemIter->next();
141 if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
142 SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
143 while ( anIter->more() ) {
144 if ( const SMDS_MeshElement* anElemNode = anIter->next() ) {
145 if ( anElemNode == aNode0 ) {
147 if ( !aNode1 ) break; // not a quadratic edge
149 else if ( anElemNode == aNode1 )
155 int aResult = std::max ( aResult0, aResult1 );
157 // TColStd_MapOfInteger aMap;
159 // SMDS_ElemIteratorPtr anIter = anEdge->nodesIterator();
160 // if ( anIter != 0 ) {
161 // while( anIter->more() ) {
162 // const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
165 // SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
166 // while( anElemIter->more() ) {
167 // const SMDS_MeshElement* anElem = anElemIter->next();
168 // if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
169 // int anId = anElem->GetID();
171 // if ( anIter->more() ) // i.e. first node
173 // else if ( aMap.Contains( anId ) )
183 gp_XYZ getNormale( const SMDS_MeshFace* theFace, bool* ok=0 )
185 int aNbNode = theFace->NbNodes();
187 gp_XYZ q1 = gpXYZ( theFace->GetNode(1)) - gpXYZ( theFace->GetNode(0));
188 gp_XYZ q2 = gpXYZ( theFace->GetNode(2)) - gpXYZ( theFace->GetNode(0));
191 gp_XYZ q3 = gpXYZ( theFace->GetNode(3)) - gpXYZ( theFace->GetNode(0));
194 double len = n.Modulus();
195 bool zeroLen = ( len <= numeric_limits<double>::min());
199 if (ok) *ok = !zeroLen;
207 using namespace SMESH::Controls;
213 //================================================================================
215 Class : NumericalFunctor
216 Description : Base class for numerical functors
218 //================================================================================
220 NumericalFunctor::NumericalFunctor():
226 void NumericalFunctor::SetMesh( const SMDS_Mesh* theMesh )
231 bool NumericalFunctor::GetPoints(const int theId,
232 TSequenceOfXYZ& theRes ) const
239 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
240 if ( !anElem || anElem->GetType() != this->GetType() )
243 return GetPoints( anElem, theRes );
246 bool NumericalFunctor::GetPoints(const SMDS_MeshElement* anElem,
247 TSequenceOfXYZ& theRes )
254 theRes.reserve( anElem->NbNodes() );
256 // Get nodes of the element
257 SMDS_ElemIteratorPtr anIter;
259 if ( anElem->IsQuadratic() ) {
260 switch ( anElem->GetType() ) {
262 anIter = dynamic_cast<const SMDS_VtkEdge*>
263 (anElem)->interlacedNodesElemIterator();
266 anIter = dynamic_cast<const SMDS_VtkFace*>
267 (anElem)->interlacedNodesElemIterator();
270 anIter = anElem->nodesIterator();
275 anIter = anElem->nodesIterator();
279 while( anIter->more() ) {
280 if ( const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>( anIter->next() ))
281 theRes.push_back( gp_XYZ( aNode->X(), aNode->Y(), aNode->Z() ) );
288 long NumericalFunctor::GetPrecision() const
293 void NumericalFunctor::SetPrecision( const long thePrecision )
295 myPrecision = thePrecision;
296 myPrecisionValue = pow( 10., (double)( myPrecision ) );
299 double NumericalFunctor::GetValue( long theId )
303 myCurrElement = myMesh->FindElement( theId );
306 if ( GetPoints( theId, P ))
307 aVal = Round( GetValue( P ));
312 double NumericalFunctor::Round( const double & aVal )
314 return ( myPrecision >= 0 ) ? floor( aVal * myPrecisionValue + 0.5 ) / myPrecisionValue : aVal;
317 //================================================================================
319 * \brief Return histogram of functor values
320 * \param nbIntervals - number of intervals
321 * \param nbEvents - number of mesh elements having values within i-th interval
322 * \param funValues - boundaries of intervals
323 * \param elements - elements to check vulue of; empty list means "of all"
324 * \param minmax - boundaries of diapason of values to divide into intervals
326 //================================================================================
328 void NumericalFunctor::GetHistogram(int nbIntervals,
329 std::vector<int>& nbEvents,
330 std::vector<double>& funValues,
331 const vector<int>& elements,
332 const double* minmax,
333 const bool isLogarithmic)
335 if ( nbIntervals < 1 ||
337 !myMesh->GetMeshInfo().NbElements( GetType() ))
339 nbEvents.resize( nbIntervals, 0 );
340 funValues.resize( nbIntervals+1 );
342 // get all values sorted
343 std::multiset< double > values;
344 if ( elements.empty() )
346 SMDS_ElemIteratorPtr elemIt = myMesh->elementsIterator(GetType());
347 while ( elemIt->more() )
348 values.insert( GetValue( elemIt->next()->GetID() ));
352 vector<int>::const_iterator id = elements.begin();
353 for ( ; id != elements.end(); ++id )
354 values.insert( GetValue( *id ));
359 funValues[0] = minmax[0];
360 funValues[nbIntervals] = minmax[1];
364 funValues[0] = *values.begin();
365 funValues[nbIntervals] = *values.rbegin();
367 // case nbIntervals == 1
368 if ( nbIntervals == 1 )
370 nbEvents[0] = values.size();
374 if (funValues.front() == funValues.back())
376 nbEvents.resize( 1 );
377 nbEvents[0] = values.size();
378 funValues[1] = funValues.back();
379 funValues.resize( 2 );
382 std::multiset< double >::iterator min = values.begin(), max;
383 for ( int i = 0; i < nbIntervals; ++i )
385 // find end value of i-th interval
386 double r = (i+1) / double(nbIntervals);
387 if (isLogarithmic && funValues.front() > 1e-07 && funValues.back() > 1e-07) {
388 double logmin = log10(funValues.front());
389 double lval = logmin + r * (log10(funValues.back()) - logmin);
390 funValues[i+1] = pow(10.0, lval);
393 funValues[i+1] = funValues.front() * (1-r) + funValues.back() * r;
396 // count values in the i-th interval if there are any
397 if ( min != values.end() && *min <= funValues[i+1] )
399 // find the first value out of the interval
400 max = values.upper_bound( funValues[i+1] ); // max is greater than funValues[i+1], or end()
401 nbEvents[i] = std::distance( min, max );
405 // add values larger than minmax[1]
406 nbEvents.back() += std::distance( min, values.end() );
409 //=======================================================================
412 Description : Functor calculating volume of a 3D element
414 //================================================================================
416 double Volume::GetValue( long theElementId )
418 if ( theElementId && myMesh ) {
419 SMDS_VolumeTool aVolumeTool;
420 if ( aVolumeTool.Set( myMesh->FindElement( theElementId )))
421 return aVolumeTool.GetSize();
426 double Volume::GetBadRate( double Value, int /*nbNodes*/ ) const
431 SMDSAbs_ElementType Volume::GetType() const
433 return SMDSAbs_Volume;
436 //=======================================================================
438 Class : MaxElementLength2D
439 Description : Functor calculating maximum length of 2D element
441 //================================================================================
443 double MaxElementLength2D::GetValue( const TSequenceOfXYZ& P )
449 if( len == 3 ) { // triangles
450 double L1 = getDistance(P( 1 ),P( 2 ));
451 double L2 = getDistance(P( 2 ),P( 3 ));
452 double L3 = getDistance(P( 3 ),P( 1 ));
453 aVal = Max(L1,Max(L2,L3));
455 else if( len == 4 ) { // quadrangles
456 double L1 = getDistance(P( 1 ),P( 2 ));
457 double L2 = getDistance(P( 2 ),P( 3 ));
458 double L3 = getDistance(P( 3 ),P( 4 ));
459 double L4 = getDistance(P( 4 ),P( 1 ));
460 double D1 = getDistance(P( 1 ),P( 3 ));
461 double D2 = getDistance(P( 2 ),P( 4 ));
462 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
464 else if( len == 6 ) { // quadratic triangles
465 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
466 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
467 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
468 aVal = Max(L1,Max(L2,L3));
470 else if( len == 8 || len == 9 ) { // quadratic quadrangles
471 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
472 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
473 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
474 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
475 double D1 = getDistance(P( 1 ),P( 5 ));
476 double D2 = getDistance(P( 3 ),P( 7 ));
477 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
480 if( myPrecision >= 0 )
482 double prec = pow( 10., (double)myPrecision );
483 aVal = floor( aVal * prec + 0.5 ) / prec;
488 double MaxElementLength2D::GetValue( long theElementId )
491 return GetPoints( theElementId, P ) ? GetValue(P) : 0.0;
494 double MaxElementLength2D::GetBadRate( double Value, int /*nbNodes*/ ) const
499 SMDSAbs_ElementType MaxElementLength2D::GetType() const
504 //=======================================================================
506 Class : MaxElementLength3D
507 Description : Functor calculating maximum length of 3D element
509 //================================================================================
511 double MaxElementLength3D::GetValue( long theElementId )
514 if( GetPoints( theElementId, P ) ) {
516 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
517 SMDSAbs_ElementType aType = aElem->GetType();
521 if( len == 4 ) { // tetras
522 double L1 = getDistance(P( 1 ),P( 2 ));
523 double L2 = getDistance(P( 2 ),P( 3 ));
524 double L3 = getDistance(P( 3 ),P( 1 ));
525 double L4 = getDistance(P( 1 ),P( 4 ));
526 double L5 = getDistance(P( 2 ),P( 4 ));
527 double L6 = getDistance(P( 3 ),P( 4 ));
528 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
531 else if( len == 5 ) { // pyramids
532 double L1 = getDistance(P( 1 ),P( 2 ));
533 double L2 = getDistance(P( 2 ),P( 3 ));
534 double L3 = getDistance(P( 3 ),P( 4 ));
535 double L4 = getDistance(P( 4 ),P( 1 ));
536 double L5 = getDistance(P( 1 ),P( 5 ));
537 double L6 = getDistance(P( 2 ),P( 5 ));
538 double L7 = getDistance(P( 3 ),P( 5 ));
539 double L8 = getDistance(P( 4 ),P( 5 ));
540 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
541 aVal = Max(aVal,Max(L7,L8));
544 else if( len == 6 ) { // pentas
545 double L1 = getDistance(P( 1 ),P( 2 ));
546 double L2 = getDistance(P( 2 ),P( 3 ));
547 double L3 = getDistance(P( 3 ),P( 1 ));
548 double L4 = getDistance(P( 4 ),P( 5 ));
549 double L5 = getDistance(P( 5 ),P( 6 ));
550 double L6 = getDistance(P( 6 ),P( 4 ));
551 double L7 = getDistance(P( 1 ),P( 4 ));
552 double L8 = getDistance(P( 2 ),P( 5 ));
553 double L9 = getDistance(P( 3 ),P( 6 ));
554 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
555 aVal = Max(aVal,Max(Max(L7,L8),L9));
558 else if( len == 8 ) { // hexas
559 double L1 = getDistance(P( 1 ),P( 2 ));
560 double L2 = getDistance(P( 2 ),P( 3 ));
561 double L3 = getDistance(P( 3 ),P( 4 ));
562 double L4 = getDistance(P( 4 ),P( 1 ));
563 double L5 = getDistance(P( 5 ),P( 6 ));
564 double L6 = getDistance(P( 6 ),P( 7 ));
565 double L7 = getDistance(P( 7 ),P( 8 ));
566 double L8 = getDistance(P( 8 ),P( 5 ));
567 double L9 = getDistance(P( 1 ),P( 5 ));
568 double L10= getDistance(P( 2 ),P( 6 ));
569 double L11= getDistance(P( 3 ),P( 7 ));
570 double L12= getDistance(P( 4 ),P( 8 ));
571 double D1 = getDistance(P( 1 ),P( 7 ));
572 double D2 = getDistance(P( 2 ),P( 8 ));
573 double D3 = getDistance(P( 3 ),P( 5 ));
574 double D4 = getDistance(P( 4 ),P( 6 ));
575 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
576 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
577 aVal = Max(aVal,Max(L11,L12));
578 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
581 else if( len == 12 ) { // hexagonal prism
582 for ( int i1 = 1; i1 < 12; ++i1 )
583 for ( int i2 = i1+1; i1 <= 12; ++i1 )
584 aVal = Max( aVal, getDistance(P( i1 ),P( i2 )));
587 else if( len == 10 ) { // quadratic tetras
588 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
589 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
590 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
591 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
592 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
593 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
594 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
597 else if( len == 13 ) { // quadratic pyramids
598 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
599 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
600 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
601 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
602 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
603 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
604 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
605 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
606 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
607 aVal = Max(aVal,Max(L7,L8));
610 else if( len == 15 ) { // quadratic pentas
611 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
612 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
613 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
614 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
615 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
616 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
617 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
618 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
619 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
620 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
621 aVal = Max(aVal,Max(Max(L7,L8),L9));
624 else if( len == 20 || len == 27 ) { // quadratic hexas
625 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
626 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
627 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
628 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
629 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
630 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
631 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
632 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
633 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
634 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
635 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
636 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
637 double D1 = getDistance(P( 1 ),P( 7 ));
638 double D2 = getDistance(P( 2 ),P( 8 ));
639 double D3 = getDistance(P( 3 ),P( 5 ));
640 double D4 = getDistance(P( 4 ),P( 6 ));
641 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
642 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
643 aVal = Max(aVal,Max(L11,L12));
644 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
647 else if( len > 1 && aElem->IsPoly() ) { // polys
648 // get the maximum distance between all pairs of nodes
649 for( int i = 1; i <= len; i++ ) {
650 for( int j = 1; j <= len; j++ ) {
651 if( j > i ) { // optimization of the loop
652 double D = getDistance( P(i), P(j) );
653 aVal = Max( aVal, D );
660 if( myPrecision >= 0 )
662 double prec = pow( 10., (double)myPrecision );
663 aVal = floor( aVal * prec + 0.5 ) / prec;
670 double MaxElementLength3D::GetBadRate( double Value, int /*nbNodes*/ ) const
675 SMDSAbs_ElementType MaxElementLength3D::GetType() const
677 return SMDSAbs_Volume;
680 //=======================================================================
683 Description : Functor for calculation of minimum angle
685 //================================================================================
687 double MinimumAngle::GetValue( const TSequenceOfXYZ& P )
694 aMin = getAngle(P( P.size() ), P( 1 ), P( 2 ));
695 aMin = Min(aMin,getAngle(P( P.size()-1 ), P( P.size() ), P( 1 )));
697 for (int i=2; i<P.size();i++){
698 double A0 = getAngle( P( i-1 ), P( i ), P( i+1 ) );
702 return aMin * 180.0 / M_PI;
705 double MinimumAngle::GetBadRate( double Value, int nbNodes ) const
707 //const double aBestAngle = PI / nbNodes;
708 const double aBestAngle = 180.0 - ( 360.0 / double(nbNodes) );
709 return ( fabs( aBestAngle - Value ));
712 SMDSAbs_ElementType MinimumAngle::GetType() const
718 //================================================================================
721 Description : Functor for calculating aspect ratio
723 //================================================================================
725 double AspectRatio::GetValue( long theId )
728 myCurrElement = myMesh->FindElement( theId );
729 if ( myCurrElement && myCurrElement->GetVtkType() == VTK_QUAD )
732 vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myCurrElement->getMeshId()]->getGrid();
733 if ( vtkCell* avtkCell = grid->GetCell( myCurrElement->getVtkId() ))
734 aVal = Round( vtkMeshQuality::QuadAspectRatio( avtkCell ));
739 if ( GetPoints( myCurrElement, P ))
740 aVal = Round( GetValue( P ));
745 double AspectRatio::GetValue( const TSequenceOfXYZ& P )
747 // According to "Mesh quality control" by Nadir Bouhamau referring to
748 // Pascal Jean Frey and Paul-Louis George. Maillages, applications aux elements finis.
749 // Hermes Science publications, Paris 1999 ISBN 2-7462-0024-4
752 int nbNodes = P.size();
757 // Compute aspect ratio
759 if ( nbNodes == 3 ) {
760 // Compute lengths of the sides
761 std::vector< double > aLen (nbNodes);
762 for ( int i = 0; i < nbNodes - 1; i++ )
763 aLen[ i ] = getDistance( P( i + 1 ), P( i + 2 ) );
764 aLen[ nbNodes - 1 ] = getDistance( P( 1 ), P( nbNodes ) );
765 // Q = alfa * h * p / S, where
767 // alfa = sqrt( 3 ) / 6
768 // h - length of the longest edge
769 // p - half perimeter
770 // S - triangle surface
771 const double alfa = sqrt( 3. ) / 6.;
772 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
773 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
774 double anArea = getArea( P( 1 ), P( 2 ), P( 3 ) );
775 if ( anArea <= theEps )
777 return alfa * maxLen * half_perimeter / anArea;
779 else if ( nbNodes == 6 ) { // quadratic triangles
780 // Compute lengths of the sides
781 std::vector< double > aLen (3);
782 aLen[0] = getDistance( P(1), P(3) );
783 aLen[1] = getDistance( P(3), P(5) );
784 aLen[2] = getDistance( P(5), P(1) );
785 // Q = alfa * h * p / S, where
787 // alfa = sqrt( 3 ) / 6
788 // h - length of the longest edge
789 // p - half perimeter
790 // S - triangle surface
791 const double alfa = sqrt( 3. ) / 6.;
792 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
793 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
794 double anArea = getArea( P(1), P(3), P(5) );
795 if ( anArea <= theEps )
797 return alfa * maxLen * half_perimeter / anArea;
799 else if( nbNodes == 4 ) { // quadrangle
800 // Compute lengths of the sides
801 std::vector< double > aLen (4);
802 aLen[0] = getDistance( P(1), P(2) );
803 aLen[1] = getDistance( P(2), P(3) );
804 aLen[2] = getDistance( P(3), P(4) );
805 aLen[3] = getDistance( P(4), P(1) );
806 // Compute lengths of the diagonals
807 std::vector< double > aDia (2);
808 aDia[0] = getDistance( P(1), P(3) );
809 aDia[1] = getDistance( P(2), P(4) );
810 // Compute areas of all triangles which can be built
811 // taking three nodes of the quadrangle
812 std::vector< double > anArea (4);
813 anArea[0] = getArea( P(1), P(2), P(3) );
814 anArea[1] = getArea( P(1), P(2), P(4) );
815 anArea[2] = getArea( P(1), P(3), P(4) );
816 anArea[3] = getArea( P(2), P(3), P(4) );
817 // Q = alpha * L * C1 / C2, where
819 // alpha = sqrt( 1/32 )
820 // L = max( L1, L2, L3, L4, D1, D2 )
821 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
822 // C2 = min( S1, S2, S3, S4 )
823 // Li - lengths of the edges
824 // Di - lengths of the diagonals
825 // Si - areas of the triangles
826 const double alpha = sqrt( 1 / 32. );
827 double L = Max( aLen[ 0 ],
831 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
832 double C1 = sqrt( ( aLen[0] * aLen[0] +
835 aLen[3] * aLen[3] ) / 4. );
836 double C2 = Min( anArea[ 0 ],
838 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
841 return alpha * L * C1 / C2;
843 else if( nbNodes == 8 || nbNodes == 9 ) { // nbNodes==8 - quadratic quadrangle
844 // Compute lengths of the sides
845 std::vector< double > aLen (4);
846 aLen[0] = getDistance( P(1), P(3) );
847 aLen[1] = getDistance( P(3), P(5) );
848 aLen[2] = getDistance( P(5), P(7) );
849 aLen[3] = getDistance( P(7), P(1) );
850 // Compute lengths of the diagonals
851 std::vector< double > aDia (2);
852 aDia[0] = getDistance( P(1), P(5) );
853 aDia[1] = getDistance( P(3), P(7) );
854 // Compute areas of all triangles which can be built
855 // taking three nodes of the quadrangle
856 std::vector< double > anArea (4);
857 anArea[0] = getArea( P(1), P(3), P(5) );
858 anArea[1] = getArea( P(1), P(3), P(7) );
859 anArea[2] = getArea( P(1), P(5), P(7) );
860 anArea[3] = getArea( P(3), P(5), P(7) );
861 // Q = alpha * L * C1 / C2, where
863 // alpha = sqrt( 1/32 )
864 // L = max( L1, L2, L3, L4, D1, D2 )
865 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
866 // C2 = min( S1, S2, S3, S4 )
867 // Li - lengths of the edges
868 // Di - lengths of the diagonals
869 // Si - areas of the triangles
870 const double alpha = sqrt( 1 / 32. );
871 double L = Max( aLen[ 0 ],
875 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
876 double C1 = sqrt( ( aLen[0] * aLen[0] +
879 aLen[3] * aLen[3] ) / 4. );
880 double C2 = Min( anArea[ 0 ],
882 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
885 return alpha * L * C1 / C2;
890 double AspectRatio::GetBadRate( double Value, int /*nbNodes*/ ) const
892 // the aspect ratio is in the range [1.0,infinity]
893 // < 1.0 = very bad, zero area
896 return ( Value < 0.9 ) ? 1000 : Value / 1000.;
899 SMDSAbs_ElementType AspectRatio::GetType() const
905 //================================================================================
907 Class : AspectRatio3D
908 Description : Functor for calculating aspect ratio
910 //================================================================================
914 inline double getHalfPerimeter(double theTria[3]){
915 return (theTria[0] + theTria[1] + theTria[2])/2.0;
918 inline double getArea(double theHalfPerim, double theTria[3]){
919 return sqrt(theHalfPerim*
920 (theHalfPerim-theTria[0])*
921 (theHalfPerim-theTria[1])*
922 (theHalfPerim-theTria[2]));
925 inline double getVolume(double theLen[6]){
926 double a2 = theLen[0]*theLen[0];
927 double b2 = theLen[1]*theLen[1];
928 double c2 = theLen[2]*theLen[2];
929 double d2 = theLen[3]*theLen[3];
930 double e2 = theLen[4]*theLen[4];
931 double f2 = theLen[5]*theLen[5];
932 double P = 4.0*a2*b2*d2;
933 double Q = a2*(b2+d2-e2)-b2*(a2+d2-f2)-d2*(a2+b2-c2);
934 double R = (b2+d2-e2)*(a2+d2-f2)*(a2+d2-f2);
935 return sqrt(P-Q+R)/12.0;
938 inline double getVolume2(double theLen[6]){
939 double a2 = theLen[0]*theLen[0];
940 double b2 = theLen[1]*theLen[1];
941 double c2 = theLen[2]*theLen[2];
942 double d2 = theLen[3]*theLen[3];
943 double e2 = theLen[4]*theLen[4];
944 double f2 = theLen[5]*theLen[5];
946 double P = a2*e2*(b2+c2+d2+f2-a2-e2);
947 double Q = b2*f2*(a2+c2+d2+e2-b2-f2);
948 double R = c2*d2*(a2+b2+e2+f2-c2-d2);
949 double S = a2*b2*d2+b2*c2*e2+a2*c2*f2+d2*e2*f2;
951 return sqrt(P+Q+R-S)/12.0;
954 inline double getVolume(const TSequenceOfXYZ& P){
955 gp_Vec aVec1( P( 2 ) - P( 1 ) );
956 gp_Vec aVec2( P( 3 ) - P( 1 ) );
957 gp_Vec aVec3( P( 4 ) - P( 1 ) );
958 gp_Vec anAreaVec( aVec1 ^ aVec2 );
959 return fabs(aVec3 * anAreaVec) / 6.0;
962 inline double getMaxHeight(double theLen[6])
964 double aHeight = std::max(theLen[0],theLen[1]);
965 aHeight = std::max(aHeight,theLen[2]);
966 aHeight = std::max(aHeight,theLen[3]);
967 aHeight = std::max(aHeight,theLen[4]);
968 aHeight = std::max(aHeight,theLen[5]);
974 double AspectRatio3D::GetValue( long theId )
977 myCurrElement = myMesh->FindElement( theId );
978 if ( myCurrElement && myCurrElement->GetVtkType() == VTK_TETRA )
980 // Action from CoTech | ACTION 31.3:
981 // EURIWARE BO: Homogenize the formulas used to calculate the Controls in SMESH to fit with
982 // those of ParaView. The library used by ParaView for those calculations can be reused in SMESH.
983 vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myCurrElement->getMeshId()]->getGrid();
984 if ( vtkCell* avtkCell = grid->GetCell( myCurrElement->getVtkId() ))
985 aVal = Round( vtkMeshQuality::TetAspectRatio( avtkCell ));
990 if ( GetPoints( myCurrElement, P ))
991 aVal = Round( GetValue( P ));
996 double AspectRatio3D::GetValue( const TSequenceOfXYZ& P )
998 double aQuality = 0.0;
999 if(myCurrElement->IsPoly()) return aQuality;
1001 int nbNodes = P.size();
1003 if(myCurrElement->IsQuadratic()) {
1004 if(nbNodes==10) nbNodes=4; // quadratic tetrahedron
1005 else if(nbNodes==13) nbNodes=5; // quadratic pyramid
1006 else if(nbNodes==15) nbNodes=6; // quadratic pentahedron
1007 else if(nbNodes==20) nbNodes=8; // quadratic hexahedron
1008 else if(nbNodes==27) nbNodes=8; // quadratic hexahedron
1009 else return aQuality;
1015 getDistance(P( 1 ),P( 2 )), // a
1016 getDistance(P( 2 ),P( 3 )), // b
1017 getDistance(P( 3 ),P( 1 )), // c
1018 getDistance(P( 2 ),P( 4 )), // d
1019 getDistance(P( 3 ),P( 4 )), // e
1020 getDistance(P( 1 ),P( 4 )) // f
1022 double aTria[4][3] = {
1023 {aLen[0],aLen[1],aLen[2]}, // abc
1024 {aLen[0],aLen[3],aLen[5]}, // adf
1025 {aLen[1],aLen[3],aLen[4]}, // bde
1026 {aLen[2],aLen[4],aLen[5]} // cef
1028 double aSumArea = 0.0;
1029 double aHalfPerimeter = getHalfPerimeter(aTria[0]);
1030 double anArea = getArea(aHalfPerimeter,aTria[0]);
1032 aHalfPerimeter = getHalfPerimeter(aTria[1]);
1033 anArea = getArea(aHalfPerimeter,aTria[1]);
1035 aHalfPerimeter = getHalfPerimeter(aTria[2]);
1036 anArea = getArea(aHalfPerimeter,aTria[2]);
1038 aHalfPerimeter = getHalfPerimeter(aTria[3]);
1039 anArea = getArea(aHalfPerimeter,aTria[3]);
1041 double aVolume = getVolume(P);
1042 //double aVolume = getVolume(aLen);
1043 double aHeight = getMaxHeight(aLen);
1044 static double aCoeff = sqrt(2.0)/12.0;
1045 if ( aVolume > DBL_MIN )
1046 aQuality = aCoeff*aHeight*aSumArea/aVolume;
1051 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 3 ),P( 5 )};
1052 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1055 gp_XYZ aXYZ[4] = {P( 1 ),P( 3 ),P( 4 ),P( 5 )};
1056 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1059 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 5 )};
1060 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1063 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 4 ),P( 5 )};
1064 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1070 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 6 )};
1071 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1074 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 3 )};
1075 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1078 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 6 )};
1079 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1082 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1083 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1086 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 6 )};
1087 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1090 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 3 )};
1091 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1097 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1098 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1101 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 4 )};
1102 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1105 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 7 )};
1106 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1109 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 8 )};
1110 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1113 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 3 )};
1114 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1117 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 4 )};
1118 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1121 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 7 )};
1122 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1125 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 8 )};
1126 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1129 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 3 )};
1130 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1133 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 4 )};
1134 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1137 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 7 )};
1138 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1141 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 8 )};
1142 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1145 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 1 )};
1146 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1149 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 2 )};
1150 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1153 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 5 )};
1154 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1157 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 6 )};
1158 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1161 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 1 )};
1162 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1165 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 2 )};
1166 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1169 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 5 )};
1170 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1173 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 6 )};
1174 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1177 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 1 )};
1178 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1181 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1182 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1185 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 5 )};
1186 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1189 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 6 )};
1190 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1193 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1194 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1197 gp_XYZ aXYZ[4] = {P( 4 ),P( 5 ),P( 8 ),P( 2 )};
1198 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1201 gp_XYZ aXYZ[4] = {P( 1 ),P( 4 ),P( 5 ),P( 3 )};
1202 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1205 gp_XYZ aXYZ[4] = {P( 3 ),P( 6 ),P( 7 ),P( 1 )};
1206 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1209 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 6 ),P( 4 )};
1210 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1213 gp_XYZ aXYZ[4] = {P( 5 ),P( 6 ),P( 8 ),P( 3 )};
1214 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1217 gp_XYZ aXYZ[4] = {P( 7 ),P( 8 ),P( 6 ),P( 1 )};
1218 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1221 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 7 )};
1222 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1225 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 2 ),P( 5 )};
1226 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1232 gp_XYZ aXYZ[8] = {P( 1 ),P( 2 ),P( 4 ),P( 5 ),P( 7 ),P( 8 ),P( 10 ),P( 11 )};
1233 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1236 gp_XYZ aXYZ[8] = {P( 2 ),P( 3 ),P( 5 ),P( 6 ),P( 8 ),P( 9 ),P( 11 ),P( 12 )};
1237 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1240 gp_XYZ aXYZ[8] = {P( 3 ),P( 4 ),P( 6 ),P( 1 ),P( 9 ),P( 10 ),P( 12 ),P( 7 )};
1241 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1244 } // switch(nbNodes)
1246 if ( nbNodes > 4 ) {
1247 // avaluate aspect ratio of quadranle faces
1248 AspectRatio aspect2D;
1249 SMDS_VolumeTool::VolumeType type = SMDS_VolumeTool::GetType( nbNodes );
1250 int nbFaces = SMDS_VolumeTool::NbFaces( type );
1251 TSequenceOfXYZ points(4);
1252 for ( int i = 0; i < nbFaces; ++i ) { // loop on faces of a volume
1253 if ( SMDS_VolumeTool::NbFaceNodes( type, i ) != 4 )
1255 const int* pInd = SMDS_VolumeTool::GetFaceNodesIndices( type, i, true );
1256 for ( int p = 0; p < 4; ++p ) // loop on nodes of a quadranle face
1257 points( p + 1 ) = P( pInd[ p ] + 1 );
1258 aQuality = std::max( aQuality, aspect2D.GetValue( points ));
1264 double AspectRatio3D::GetBadRate( double Value, int /*nbNodes*/ ) const
1266 // the aspect ratio is in the range [1.0,infinity]
1269 return Value / 1000.;
1272 SMDSAbs_ElementType AspectRatio3D::GetType() const
1274 return SMDSAbs_Volume;
1278 //================================================================================
1281 Description : Functor for calculating warping
1283 //================================================================================
1285 double Warping::GetValue( const TSequenceOfXYZ& P )
1287 if ( P.size() != 4 )
1290 gp_XYZ G = ( P( 1 ) + P( 2 ) + P( 3 ) + P( 4 ) ) / 4.;
1292 double A1 = ComputeA( P( 1 ), P( 2 ), P( 3 ), G );
1293 double A2 = ComputeA( P( 2 ), P( 3 ), P( 4 ), G );
1294 double A3 = ComputeA( P( 3 ), P( 4 ), P( 1 ), G );
1295 double A4 = ComputeA( P( 4 ), P( 1 ), P( 2 ), G );
1297 return Max( Max( A1, A2 ), Max( A3, A4 ) );
1300 double Warping::ComputeA( const gp_XYZ& thePnt1,
1301 const gp_XYZ& thePnt2,
1302 const gp_XYZ& thePnt3,
1303 const gp_XYZ& theG ) const
1305 double aLen1 = gp_Pnt( thePnt1 ).Distance( gp_Pnt( thePnt2 ) );
1306 double aLen2 = gp_Pnt( thePnt2 ).Distance( gp_Pnt( thePnt3 ) );
1307 double L = Min( aLen1, aLen2 ) * 0.5;
1311 gp_XYZ GI = ( thePnt2 + thePnt1 ) / 2. - theG;
1312 gp_XYZ GJ = ( thePnt3 + thePnt2 ) / 2. - theG;
1313 gp_XYZ N = GI.Crossed( GJ );
1315 if ( N.Modulus() < gp::Resolution() )
1320 double H = ( thePnt2 - theG ).Dot( N );
1321 return asin( fabs( H / L ) ) * 180. / M_PI;
1324 double Warping::GetBadRate( double Value, int /*nbNodes*/ ) const
1326 // the warp is in the range [0.0,PI/2]
1327 // 0.0 = good (no warp)
1328 // PI/2 = bad (face pliee)
1332 SMDSAbs_ElementType Warping::GetType() const
1334 return SMDSAbs_Face;
1338 //================================================================================
1341 Description : Functor for calculating taper
1343 //================================================================================
1345 double Taper::GetValue( const TSequenceOfXYZ& P )
1347 if ( P.size() != 4 )
1351 double J1 = getArea( P( 4 ), P( 1 ), P( 2 ) ) / 2.;
1352 double J2 = getArea( P( 3 ), P( 1 ), P( 2 ) ) / 2.;
1353 double J3 = getArea( P( 2 ), P( 3 ), P( 4 ) ) / 2.;
1354 double J4 = getArea( P( 3 ), P( 4 ), P( 1 ) ) / 2.;
1356 double JA = 0.25 * ( J1 + J2 + J3 + J4 );
1360 double T1 = fabs( ( J1 - JA ) / JA );
1361 double T2 = fabs( ( J2 - JA ) / JA );
1362 double T3 = fabs( ( J3 - JA ) / JA );
1363 double T4 = fabs( ( J4 - JA ) / JA );
1365 return Max( Max( T1, T2 ), Max( T3, T4 ) );
1368 double Taper::GetBadRate( double Value, int /*nbNodes*/ ) const
1370 // the taper is in the range [0.0,1.0]
1371 // 0.0 = good (no taper)
1372 // 1.0 = bad (les cotes opposes sont allignes)
1376 SMDSAbs_ElementType Taper::GetType() const
1378 return SMDSAbs_Face;
1381 //================================================================================
1384 Description : Functor for calculating skew in degrees
1386 //================================================================================
1388 static inline double skewAngle( const gp_XYZ& p1, const gp_XYZ& p2, const gp_XYZ& p3 )
1390 gp_XYZ p12 = ( p2 + p1 ) / 2.;
1391 gp_XYZ p23 = ( p3 + p2 ) / 2.;
1392 gp_XYZ p31 = ( p3 + p1 ) / 2.;
1394 gp_Vec v1( p31 - p2 ), v2( p12 - p23 );
1396 return v1.Magnitude() < gp::Resolution() || v2.Magnitude() < gp::Resolution() ? 0. : v1.Angle( v2 );
1399 double Skew::GetValue( const TSequenceOfXYZ& P )
1401 if ( P.size() != 3 && P.size() != 4 )
1405 static double PI2 = M_PI / 2.;
1406 if ( P.size() == 3 )
1408 double A0 = fabs( PI2 - skewAngle( P( 3 ), P( 1 ), P( 2 ) ) );
1409 double A1 = fabs( PI2 - skewAngle( P( 1 ), P( 2 ), P( 3 ) ) );
1410 double A2 = fabs( PI2 - skewAngle( P( 2 ), P( 3 ), P( 1 ) ) );
1412 return Max( A0, Max( A1, A2 ) ) * 180. / M_PI;
1416 gp_XYZ p12 = ( P( 1 ) + P( 2 ) ) / 2.;
1417 gp_XYZ p23 = ( P( 2 ) + P( 3 ) ) / 2.;
1418 gp_XYZ p34 = ( P( 3 ) + P( 4 ) ) / 2.;
1419 gp_XYZ p41 = ( P( 4 ) + P( 1 ) ) / 2.;
1421 gp_Vec v1( p34 - p12 ), v2( p23 - p41 );
1422 double A = v1.Magnitude() <= gp::Resolution() || v2.Magnitude() <= gp::Resolution()
1423 ? 0. : fabs( PI2 - v1.Angle( v2 ) );
1429 return A * 180. / M_PI;
1433 double Skew::GetBadRate( double Value, int /*nbNodes*/ ) const
1435 // the skew is in the range [0.0,PI/2].
1441 SMDSAbs_ElementType Skew::GetType() const
1443 return SMDSAbs_Face;
1447 //================================================================================
1450 Description : Functor for calculating area
1452 //================================================================================
1454 double Area::GetValue( const TSequenceOfXYZ& P )
1457 if ( P.size() > 2 ) {
1458 gp_Vec aVec1( P(2) - P(1) );
1459 gp_Vec aVec2( P(3) - P(1) );
1460 gp_Vec SumVec = aVec1 ^ aVec2;
1461 for (int i=4; i<=P.size(); i++) {
1462 gp_Vec aVec1( P(i-1) - P(1) );
1463 gp_Vec aVec2( P(i) - P(1) );
1464 gp_Vec tmp = aVec1 ^ aVec2;
1467 val = SumVec.Magnitude() * 0.5;
1472 double Area::GetBadRate( double Value, int /*nbNodes*/ ) const
1474 // meaningless as it is not a quality control functor
1478 SMDSAbs_ElementType Area::GetType() const
1480 return SMDSAbs_Face;
1483 //================================================================================
1486 Description : Functor for calculating length of edge
1488 //================================================================================
1490 double Length::GetValue( const TSequenceOfXYZ& P )
1492 switch ( P.size() ) {
1493 case 2: return getDistance( P( 1 ), P( 2 ) );
1494 case 3: return getDistance( P( 1 ), P( 2 ) ) + getDistance( P( 2 ), P( 3 ) );
1499 double Length::GetBadRate( double Value, int /*nbNodes*/ ) const
1501 // meaningless as it is not quality control functor
1505 SMDSAbs_ElementType Length::GetType() const
1507 return SMDSAbs_Edge;
1510 //================================================================================
1513 Description : Functor for calculating length of edge
1515 //================================================================================
1517 double Length2D::GetValue( long theElementId)
1521 //cout<<"Length2D::GetValue"<<endl;
1522 if (GetPoints(theElementId,P)){
1523 //for(int jj=1; jj<=P.size(); jj++)
1524 // cout<<"jj="<<jj<<" P("<<P(jj).X()<<","<<P(jj).Y()<<","<<P(jj).Z()<<")"<<endl;
1526 double aVal;// = GetValue( P );
1527 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
1528 SMDSAbs_ElementType aType = aElem->GetType();
1537 aVal = getDistance( P( 1 ), P( 2 ) );
1540 else if (len == 3){ // quadratic edge
1541 aVal = getDistance(P( 1 ),P( 3 )) + getDistance(P( 3 ),P( 2 ));
1545 if (len == 3){ // triangles
1546 double L1 = getDistance(P( 1 ),P( 2 ));
1547 double L2 = getDistance(P( 2 ),P( 3 ));
1548 double L3 = getDistance(P( 3 ),P( 1 ));
1549 aVal = Max(L1,Max(L2,L3));
1552 else if (len == 4){ // quadrangles
1553 double L1 = getDistance(P( 1 ),P( 2 ));
1554 double L2 = getDistance(P( 2 ),P( 3 ));
1555 double L3 = getDistance(P( 3 ),P( 4 ));
1556 double L4 = getDistance(P( 4 ),P( 1 ));
1557 aVal = Max(Max(L1,L2),Max(L3,L4));
1560 if (len == 6){ // quadratic triangles
1561 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1562 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1563 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
1564 aVal = Max(L1,Max(L2,L3));
1565 //cout<<"L1="<<L1<<" L2="<<L2<<"L3="<<L3<<" aVal="<<aVal<<endl;
1568 else if (len == 8){ // quadratic quadrangles
1569 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1570 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1571 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
1572 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
1573 aVal = Max(Max(L1,L2),Max(L3,L4));
1576 case SMDSAbs_Volume:
1577 if (len == 4){ // tetraidrs
1578 double L1 = getDistance(P( 1 ),P( 2 ));
1579 double L2 = getDistance(P( 2 ),P( 3 ));
1580 double L3 = getDistance(P( 3 ),P( 1 ));
1581 double L4 = getDistance(P( 1 ),P( 4 ));
1582 double L5 = getDistance(P( 2 ),P( 4 ));
1583 double L6 = getDistance(P( 3 ),P( 4 ));
1584 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1587 else if (len == 5){ // piramids
1588 double L1 = getDistance(P( 1 ),P( 2 ));
1589 double L2 = getDistance(P( 2 ),P( 3 ));
1590 double L3 = getDistance(P( 3 ),P( 4 ));
1591 double L4 = getDistance(P( 4 ),P( 1 ));
1592 double L5 = getDistance(P( 1 ),P( 5 ));
1593 double L6 = getDistance(P( 2 ),P( 5 ));
1594 double L7 = getDistance(P( 3 ),P( 5 ));
1595 double L8 = getDistance(P( 4 ),P( 5 ));
1597 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1598 aVal = Max(aVal,Max(L7,L8));
1601 else if (len == 6){ // pentaidres
1602 double L1 = getDistance(P( 1 ),P( 2 ));
1603 double L2 = getDistance(P( 2 ),P( 3 ));
1604 double L3 = getDistance(P( 3 ),P( 1 ));
1605 double L4 = getDistance(P( 4 ),P( 5 ));
1606 double L5 = getDistance(P( 5 ),P( 6 ));
1607 double L6 = getDistance(P( 6 ),P( 4 ));
1608 double L7 = getDistance(P( 1 ),P( 4 ));
1609 double L8 = getDistance(P( 2 ),P( 5 ));
1610 double L9 = getDistance(P( 3 ),P( 6 ));
1612 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1613 aVal = Max(aVal,Max(Max(L7,L8),L9));
1616 else if (len == 8){ // hexaider
1617 double L1 = getDistance(P( 1 ),P( 2 ));
1618 double L2 = getDistance(P( 2 ),P( 3 ));
1619 double L3 = getDistance(P( 3 ),P( 4 ));
1620 double L4 = getDistance(P( 4 ),P( 1 ));
1621 double L5 = getDistance(P( 5 ),P( 6 ));
1622 double L6 = getDistance(P( 6 ),P( 7 ));
1623 double L7 = getDistance(P( 7 ),P( 8 ));
1624 double L8 = getDistance(P( 8 ),P( 5 ));
1625 double L9 = getDistance(P( 1 ),P( 5 ));
1626 double L10= getDistance(P( 2 ),P( 6 ));
1627 double L11= getDistance(P( 3 ),P( 7 ));
1628 double L12= getDistance(P( 4 ),P( 8 ));
1630 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1631 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1632 aVal = Max(aVal,Max(L11,L12));
1637 if (len == 10){ // quadratic tetraidrs
1638 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
1639 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
1640 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
1641 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1642 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
1643 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
1644 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1647 else if (len == 13){ // quadratic piramids
1648 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
1649 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
1650 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1651 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1652 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1653 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
1654 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
1655 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
1656 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1657 aVal = Max(aVal,Max(L7,L8));
1660 else if (len == 15){ // quadratic pentaidres
1661 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
1662 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
1663 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1664 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1665 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
1666 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
1667 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
1668 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
1669 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
1670 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1671 aVal = Max(aVal,Max(Max(L7,L8),L9));
1674 else if (len == 20){ // quadratic hexaider
1675 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
1676 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
1677 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
1678 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
1679 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
1680 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
1681 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
1682 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
1683 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
1684 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
1685 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
1686 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
1687 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1688 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1689 aVal = Max(aVal,Max(L11,L12));
1701 if ( myPrecision >= 0 )
1703 double prec = pow( 10., (double)( myPrecision ) );
1704 aVal = floor( aVal * prec + 0.5 ) / prec;
1713 double Length2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1715 // meaningless as it is not quality control functor
1719 SMDSAbs_ElementType Length2D::GetType() const
1721 return SMDSAbs_Face;
1724 Length2D::Value::Value(double theLength,long thePntId1, long thePntId2):
1727 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1728 if(thePntId1 > thePntId2){
1729 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1733 bool Length2D::Value::operator<(const Length2D::Value& x) const{
1734 if(myPntId[0] < x.myPntId[0]) return true;
1735 if(myPntId[0] == x.myPntId[0])
1736 if(myPntId[1] < x.myPntId[1]) return true;
1740 void Length2D::GetValues(TValues& theValues){
1742 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1743 for(; anIter->more(); ){
1744 const SMDS_MeshFace* anElem = anIter->next();
1746 if(anElem->IsQuadratic()) {
1747 const SMDS_VtkFace* F =
1748 dynamic_cast<const SMDS_VtkFace*>(anElem);
1749 // use special nodes iterator
1750 SMDS_ElemIteratorPtr anIter = F->interlacedNodesElemIterator();
1755 const SMDS_MeshElement* aNode;
1757 aNode = anIter->next();
1758 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1759 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1760 aNodeId[0] = aNodeId[1] = aNode->GetID();
1763 for(; anIter->more(); ){
1764 const SMDS_MeshNode* N1 = static_cast<const SMDS_MeshNode*> (anIter->next());
1765 P[2] = gp_Pnt(N1->X(),N1->Y(),N1->Z());
1766 aNodeId[2] = N1->GetID();
1767 aLength = P[1].Distance(P[2]);
1768 if(!anIter->more()) break;
1769 const SMDS_MeshNode* N2 = static_cast<const SMDS_MeshNode*> (anIter->next());
1770 P[3] = gp_Pnt(N2->X(),N2->Y(),N2->Z());
1771 aNodeId[3] = N2->GetID();
1772 aLength += P[2].Distance(P[3]);
1773 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1774 Value aValue2(aLength,aNodeId[2],aNodeId[3]);
1776 aNodeId[1] = aNodeId[3];
1777 theValues.insert(aValue1);
1778 theValues.insert(aValue2);
1780 aLength += P[2].Distance(P[0]);
1781 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1782 Value aValue2(aLength,aNodeId[2],aNodeId[0]);
1783 theValues.insert(aValue1);
1784 theValues.insert(aValue2);
1787 SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
1792 const SMDS_MeshElement* aNode;
1793 if(aNodesIter->more()){
1794 aNode = aNodesIter->next();
1795 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1796 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1797 aNodeId[0] = aNodeId[1] = aNode->GetID();
1800 for(; aNodesIter->more(); ){
1801 aNode = aNodesIter->next();
1802 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1803 long anId = aNode->GetID();
1805 P[2] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1807 aLength = P[1].Distance(P[2]);
1809 Value aValue(aLength,aNodeId[1],anId);
1812 theValues.insert(aValue);
1815 aLength = P[0].Distance(P[1]);
1817 Value aValue(aLength,aNodeId[0],aNodeId[1]);
1818 theValues.insert(aValue);
1823 //================================================================================
1825 Class : MultiConnection
1826 Description : Functor for calculating number of faces conneted to the edge
1828 //================================================================================
1830 double MultiConnection::GetValue( const TSequenceOfXYZ& P )
1834 double MultiConnection::GetValue( long theId )
1836 return getNbMultiConnection( myMesh, theId );
1839 double MultiConnection::GetBadRate( double Value, int /*nbNodes*/ ) const
1841 // meaningless as it is not quality control functor
1845 SMDSAbs_ElementType MultiConnection::GetType() const
1847 return SMDSAbs_Edge;
1850 //================================================================================
1852 Class : MultiConnection2D
1853 Description : Functor for calculating number of faces conneted to the edge
1855 //================================================================================
1857 double MultiConnection2D::GetValue( const TSequenceOfXYZ& P )
1862 double MultiConnection2D::GetValue( long theElementId )
1866 const SMDS_MeshElement* aFaceElem = myMesh->FindElement(theElementId);
1867 SMDSAbs_ElementType aType = aFaceElem->GetType();
1872 int i = 0, len = aFaceElem->NbNodes();
1873 SMDS_ElemIteratorPtr anIter = aFaceElem->nodesIterator();
1876 const SMDS_MeshNode *aNode, *aNode0;
1877 TColStd_MapOfInteger aMap, aMapPrev;
1879 for (i = 0; i <= len; i++) {
1884 if (anIter->more()) {
1885 aNode = (SMDS_MeshNode*)anIter->next();
1893 if (i == 0) aNode0 = aNode;
1895 SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
1896 while (anElemIter->more()) {
1897 const SMDS_MeshElement* anElem = anElemIter->next();
1898 if (anElem != 0 && anElem->GetType() == SMDSAbs_Face) {
1899 int anId = anElem->GetID();
1902 if (aMapPrev.Contains(anId)) {
1907 aResult = Max(aResult, aNb);
1918 double MultiConnection2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1920 // meaningless as it is not quality control functor
1924 SMDSAbs_ElementType MultiConnection2D::GetType() const
1926 return SMDSAbs_Face;
1929 MultiConnection2D::Value::Value(long thePntId1, long thePntId2)
1931 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1932 if(thePntId1 > thePntId2){
1933 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1937 bool MultiConnection2D::Value::operator<(const MultiConnection2D::Value& x) const{
1938 if(myPntId[0] < x.myPntId[0]) return true;
1939 if(myPntId[0] == x.myPntId[0])
1940 if(myPntId[1] < x.myPntId[1]) return true;
1944 void MultiConnection2D::GetValues(MValues& theValues){
1945 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1946 for(; anIter->more(); ){
1947 const SMDS_MeshFace* anElem = anIter->next();
1948 SMDS_ElemIteratorPtr aNodesIter;
1949 if ( anElem->IsQuadratic() )
1950 aNodesIter = dynamic_cast<const SMDS_VtkFace*>
1951 (anElem)->interlacedNodesElemIterator();
1953 aNodesIter = anElem->nodesIterator();
1956 //int aNbConnects=0;
1957 const SMDS_MeshNode* aNode0;
1958 const SMDS_MeshNode* aNode1;
1959 const SMDS_MeshNode* aNode2;
1960 if(aNodesIter->more()){
1961 aNode0 = (SMDS_MeshNode*) aNodesIter->next();
1963 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode1;
1964 aNodeId[0] = aNodeId[1] = aNodes->GetID();
1966 for(; aNodesIter->more(); ) {
1967 aNode2 = (SMDS_MeshNode*) aNodesIter->next();
1968 long anId = aNode2->GetID();
1971 Value aValue(aNodeId[1],aNodeId[2]);
1972 MValues::iterator aItr = theValues.find(aValue);
1973 if (aItr != theValues.end()){
1978 theValues[aValue] = 1;
1981 //cout << "NodeIds: "<<aNodeId[1]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1982 aNodeId[1] = aNodeId[2];
1985 Value aValue(aNodeId[0],aNodeId[2]);
1986 MValues::iterator aItr = theValues.find(aValue);
1987 if (aItr != theValues.end()) {
1992 theValues[aValue] = 1;
1995 //cout << "NodeIds: "<<aNodeId[0]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
2000 //================================================================================
2002 Class : BallDiameter
2003 Description : Functor returning diameter of a ball element
2005 //================================================================================
2007 double BallDiameter::GetValue( long theId )
2009 double diameter = 0;
2011 if ( const SMDS_BallElement* ball =
2012 dynamic_cast<const SMDS_BallElement*>( myMesh->FindElement( theId )))
2014 diameter = ball->GetDiameter();
2019 double BallDiameter::GetBadRate( double Value, int /*nbNodes*/ ) const
2021 // meaningless as it is not a quality control functor
2025 SMDSAbs_ElementType BallDiameter::GetType() const
2027 return SMDSAbs_Ball;
2035 //================================================================================
2037 Class : BadOrientedVolume
2038 Description : Predicate bad oriented volumes
2040 //================================================================================
2042 BadOrientedVolume::BadOrientedVolume()
2047 void BadOrientedVolume::SetMesh( const SMDS_Mesh* theMesh )
2052 bool BadOrientedVolume::IsSatisfy( long theId )
2057 SMDS_VolumeTool vTool( myMesh->FindElement( theId ));
2058 return !vTool.IsForward();
2061 SMDSAbs_ElementType BadOrientedVolume::GetType() const
2063 return SMDSAbs_Volume;
2067 Class : BareBorderVolume
2070 bool BareBorderVolume::IsSatisfy(long theElementId )
2072 SMDS_VolumeTool myTool;
2073 if ( myTool.Set( myMesh->FindElement(theElementId)))
2075 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2076 if ( myTool.IsFreeFace( iF ))
2078 const SMDS_MeshNode** n = myTool.GetFaceNodes(iF);
2079 vector< const SMDS_MeshNode*> nodes( n, n+myTool.NbFaceNodes(iF));
2080 if ( !myMesh->FindElement( nodes, SMDSAbs_Face, /*Nomedium=*/false))
2087 //================================================================================
2089 Class : BareBorderFace
2091 //================================================================================
2093 bool BareBorderFace::IsSatisfy(long theElementId )
2096 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2098 if ( face->GetType() == SMDSAbs_Face )
2100 int nbN = face->NbCornerNodes();
2101 for ( int i = 0; i < nbN && !ok; ++i )
2103 // check if a link is shared by another face
2104 const SMDS_MeshNode* n1 = face->GetNode( i );
2105 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2106 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2107 bool isShared = false;
2108 while ( !isShared && fIt->more() )
2110 const SMDS_MeshElement* f = fIt->next();
2111 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2115 const int iQuad = face->IsQuadratic();
2116 myLinkNodes.resize( 2 + iQuad);
2117 myLinkNodes[0] = n1;
2118 myLinkNodes[1] = n2;
2120 myLinkNodes[2] = face->GetNode( i+nbN );
2121 ok = !myMesh->FindElement( myLinkNodes, SMDSAbs_Edge, /*noMedium=*/false);
2129 //================================================================================
2131 Class : OverConstrainedVolume
2133 //================================================================================
2135 bool OverConstrainedVolume::IsSatisfy(long theElementId )
2137 // An element is over-constrained if it has N-1 free borders where
2138 // N is the number of edges/faces for a 2D/3D element.
2139 SMDS_VolumeTool myTool;
2140 if ( myTool.Set( myMesh->FindElement(theElementId)))
2142 int nbSharedFaces = 0;
2143 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2144 if ( !myTool.IsFreeFace( iF ) && ++nbSharedFaces > 1 )
2146 return ( nbSharedFaces == 1 );
2151 //================================================================================
2153 Class : OverConstrainedFace
2155 //================================================================================
2157 bool OverConstrainedFace::IsSatisfy(long theElementId )
2159 // An element is over-constrained if it has N-1 free borders where
2160 // N is the number of edges/faces for a 2D/3D element.
2161 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2162 if ( face->GetType() == SMDSAbs_Face )
2164 int nbSharedBorders = 0;
2165 int nbN = face->NbCornerNodes();
2166 for ( int i = 0; i < nbN; ++i )
2168 // check if a link is shared by another face
2169 const SMDS_MeshNode* n1 = face->GetNode( i );
2170 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2171 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2172 bool isShared = false;
2173 while ( !isShared && fIt->more() )
2175 const SMDS_MeshElement* f = fIt->next();
2176 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2178 if ( isShared && ++nbSharedBorders > 1 )
2181 return ( nbSharedBorders == 1 );
2186 //================================================================================
2188 Class : CoincidentNodes
2189 Description : Predicate of Coincident nodes
2191 //================================================================================
2193 CoincidentNodes::CoincidentNodes()
2198 bool CoincidentNodes::IsSatisfy( long theElementId )
2200 return myCoincidentIDs.Contains( theElementId );
2203 SMDSAbs_ElementType CoincidentNodes::GetType() const
2205 return SMDSAbs_Node;
2208 void CoincidentNodes::SetMesh( const SMDS_Mesh* theMesh )
2210 myMeshModifTracer.SetMesh( theMesh );
2211 if ( myMeshModifTracer.IsMeshModified() )
2213 TIDSortedNodeSet nodesToCheck;
2214 SMDS_NodeIteratorPtr nIt = theMesh->nodesIterator(/*idInceasingOrder=*/true);
2215 while ( nIt->more() )
2216 nodesToCheck.insert( nodesToCheck.end(), nIt->next() );
2218 list< list< const SMDS_MeshNode*> > nodeGroups;
2219 SMESH_OctreeNode::FindCoincidentNodes ( nodesToCheck, &nodeGroups, myToler );
2221 myCoincidentIDs.Clear();
2222 list< list< const SMDS_MeshNode*> >::iterator groupIt = nodeGroups.begin();
2223 for ( ; groupIt != nodeGroups.end(); ++groupIt )
2225 list< const SMDS_MeshNode*>& coincNodes = *groupIt;
2226 list< const SMDS_MeshNode*>::iterator n = coincNodes.begin();
2227 for ( ; n != coincNodes.end(); ++n )
2228 myCoincidentIDs.Add( (*n)->GetID() );
2233 //================================================================================
2235 Class : CoincidentElements
2236 Description : Predicate of Coincident Elements
2237 Note : This class is suitable only for visualization of Coincident Elements
2239 //================================================================================
2241 CoincidentElements::CoincidentElements()
2246 void CoincidentElements::SetMesh( const SMDS_Mesh* theMesh )
2251 bool CoincidentElements::IsSatisfy( long theElementId )
2253 if ( !myMesh ) return false;
2255 if ( const SMDS_MeshElement* e = myMesh->FindElement( theElementId ))
2257 if ( e->GetType() != GetType() ) return false;
2258 set< const SMDS_MeshNode* > elemNodes( e->begin_nodes(), e->end_nodes() );
2259 const int nbNodes = e->NbNodes();
2260 SMDS_ElemIteratorPtr invIt = (*elemNodes.begin())->GetInverseElementIterator( GetType() );
2261 while ( invIt->more() )
2263 const SMDS_MeshElement* e2 = invIt->next();
2264 if ( e2 == e || e2->NbNodes() != nbNodes ) continue;
2266 bool sameNodes = true;
2267 for ( size_t i = 0; i < elemNodes.size() && sameNodes; ++i )
2268 sameNodes = ( elemNodes.count( e2->GetNode( i )));
2276 SMDSAbs_ElementType CoincidentElements1D::GetType() const
2278 return SMDSAbs_Edge;
2280 SMDSAbs_ElementType CoincidentElements2D::GetType() const
2282 return SMDSAbs_Face;
2284 SMDSAbs_ElementType CoincidentElements3D::GetType() const
2286 return SMDSAbs_Volume;
2290 //================================================================================
2293 Description : Predicate for free borders
2295 //================================================================================
2297 FreeBorders::FreeBorders()
2302 void FreeBorders::SetMesh( const SMDS_Mesh* theMesh )
2307 bool FreeBorders::IsSatisfy( long theId )
2309 return getNbMultiConnection( myMesh, theId ) == 1;
2312 SMDSAbs_ElementType FreeBorders::GetType() const
2314 return SMDSAbs_Edge;
2318 //================================================================================
2321 Description : Predicate for free Edges
2323 //================================================================================
2325 FreeEdges::FreeEdges()
2330 void FreeEdges::SetMesh( const SMDS_Mesh* theMesh )
2335 bool FreeEdges::IsFreeEdge( const SMDS_MeshNode** theNodes, const int theFaceId )
2337 TColStd_MapOfInteger aMap;
2338 for ( int i = 0; i < 2; i++ )
2340 SMDS_ElemIteratorPtr anElemIter = theNodes[ i ]->GetInverseElementIterator(SMDSAbs_Face);
2341 while( anElemIter->more() )
2343 if ( const SMDS_MeshElement* anElem = anElemIter->next())
2345 const int anId = anElem->GetID();
2346 if ( anId != theFaceId && !aMap.Add( anId ))
2354 bool FreeEdges::IsSatisfy( long theId )
2359 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2360 if ( aFace == 0 || aFace->GetType() != SMDSAbs_Face || aFace->NbNodes() < 3 )
2363 SMDS_ElemIteratorPtr anIter;
2364 if ( aFace->IsQuadratic() ) {
2365 anIter = dynamic_cast<const SMDS_VtkFace*>
2366 (aFace)->interlacedNodesElemIterator();
2369 anIter = aFace->nodesIterator();
2374 int i = 0, nbNodes = aFace->NbNodes();
2375 std::vector <const SMDS_MeshNode*> aNodes( nbNodes+1 );
2376 while( anIter->more() )
2378 const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
2381 aNodes[ i++ ] = aNode;
2383 aNodes[ nbNodes ] = aNodes[ 0 ];
2385 for ( i = 0; i < nbNodes; i++ )
2386 if ( IsFreeEdge( &aNodes[ i ], theId ) )
2392 SMDSAbs_ElementType FreeEdges::GetType() const
2394 return SMDSAbs_Face;
2397 FreeEdges::Border::Border(long theElemId, long thePntId1, long thePntId2):
2400 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
2401 if(thePntId1 > thePntId2){
2402 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
2406 bool FreeEdges::Border::operator<(const FreeEdges::Border& x) const{
2407 if(myPntId[0] < x.myPntId[0]) return true;
2408 if(myPntId[0] == x.myPntId[0])
2409 if(myPntId[1] < x.myPntId[1]) return true;
2413 inline void UpdateBorders(const FreeEdges::Border& theBorder,
2414 FreeEdges::TBorders& theRegistry,
2415 FreeEdges::TBorders& theContainer)
2417 if(theRegistry.find(theBorder) == theRegistry.end()){
2418 theRegistry.insert(theBorder);
2419 theContainer.insert(theBorder);
2421 theContainer.erase(theBorder);
2425 void FreeEdges::GetBoreders(TBorders& theBorders)
2428 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
2429 for(; anIter->more(); ){
2430 const SMDS_MeshFace* anElem = anIter->next();
2431 long anElemId = anElem->GetID();
2432 SMDS_ElemIteratorPtr aNodesIter;
2433 if ( anElem->IsQuadratic() )
2434 aNodesIter = static_cast<const SMDS_VtkFace*>(anElem)->
2435 interlacedNodesElemIterator();
2437 aNodesIter = anElem->nodesIterator();
2439 const SMDS_MeshElement* aNode;
2440 if(aNodesIter->more()){
2441 aNode = aNodesIter->next();
2442 aNodeId[0] = aNodeId[1] = aNode->GetID();
2444 for(; aNodesIter->more(); ){
2445 aNode = aNodesIter->next();
2446 long anId = aNode->GetID();
2447 Border aBorder(anElemId,aNodeId[1],anId);
2449 UpdateBorders(aBorder,aRegistry,theBorders);
2451 Border aBorder(anElemId,aNodeId[0],aNodeId[1]);
2452 UpdateBorders(aBorder,aRegistry,theBorders);
2456 //================================================================================
2459 Description : Predicate for free nodes
2461 //================================================================================
2463 FreeNodes::FreeNodes()
2468 void FreeNodes::SetMesh( const SMDS_Mesh* theMesh )
2473 bool FreeNodes::IsSatisfy( long theNodeId )
2475 const SMDS_MeshNode* aNode = myMesh->FindNode( theNodeId );
2479 return (aNode->NbInverseElements() < 1);
2482 SMDSAbs_ElementType FreeNodes::GetType() const
2484 return SMDSAbs_Node;
2488 //================================================================================
2491 Description : Predicate for free faces
2493 //================================================================================
2495 FreeFaces::FreeFaces()
2500 void FreeFaces::SetMesh( const SMDS_Mesh* theMesh )
2505 bool FreeFaces::IsSatisfy( long theId )
2507 if (!myMesh) return false;
2508 // check that faces nodes refers to less than two common volumes
2509 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2510 if ( !aFace || aFace->GetType() != SMDSAbs_Face )
2513 int nbNode = aFace->NbNodes();
2515 // collect volumes check that number of volumss with count equal nbNode not less than 2
2516 typedef map< SMDS_MeshElement*, int > TMapOfVolume; // map of volume counters
2517 typedef map< SMDS_MeshElement*, int >::iterator TItrMapOfVolume; // iterator
2518 TMapOfVolume mapOfVol;
2520 SMDS_ElemIteratorPtr nodeItr = aFace->nodesIterator();
2521 while ( nodeItr->more() ) {
2522 const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>(nodeItr->next());
2523 if ( !aNode ) continue;
2524 SMDS_ElemIteratorPtr volItr = aNode->GetInverseElementIterator(SMDSAbs_Volume);
2525 while ( volItr->more() ) {
2526 SMDS_MeshElement* aVol = (SMDS_MeshElement*)volItr->next();
2527 TItrMapOfVolume itr = mapOfVol.insert(make_pair(aVol, 0)).first;
2532 TItrMapOfVolume volItr = mapOfVol.begin();
2533 TItrMapOfVolume volEnd = mapOfVol.end();
2534 for ( ; volItr != volEnd; ++volItr )
2535 if ( (*volItr).second >= nbNode )
2537 // face is not free if number of volumes constructed on thier nodes more than one
2541 SMDSAbs_ElementType FreeFaces::GetType() const
2543 return SMDSAbs_Face;
2546 //================================================================================
2548 Class : LinearOrQuadratic
2549 Description : Predicate to verify whether a mesh element is linear
2551 //================================================================================
2553 LinearOrQuadratic::LinearOrQuadratic()
2558 void LinearOrQuadratic::SetMesh( const SMDS_Mesh* theMesh )
2563 bool LinearOrQuadratic::IsSatisfy( long theId )
2565 if (!myMesh) return false;
2566 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2567 if ( !anElem || (myType != SMDSAbs_All && anElem->GetType() != myType) )
2569 return (!anElem->IsQuadratic());
2572 void LinearOrQuadratic::SetType( SMDSAbs_ElementType theType )
2577 SMDSAbs_ElementType LinearOrQuadratic::GetType() const
2582 //================================================================================
2585 Description : Functor for check color of group to whic mesh element belongs to
2587 //================================================================================
2589 GroupColor::GroupColor()
2593 bool GroupColor::IsSatisfy( long theId )
2595 return (myIDs.find( theId ) != myIDs.end());
2598 void GroupColor::SetType( SMDSAbs_ElementType theType )
2603 SMDSAbs_ElementType GroupColor::GetType() const
2608 static bool isEqual( const Quantity_Color& theColor1,
2609 const Quantity_Color& theColor2 )
2611 // tolerance to compare colors
2612 const double tol = 5*1e-3;
2613 return ( fabs( theColor1.Red() - theColor2.Red() ) < tol &&
2614 fabs( theColor1.Green() - theColor2.Green() ) < tol &&
2615 fabs( theColor1.Blue() - theColor2.Blue() ) < tol );
2619 void GroupColor::SetMesh( const SMDS_Mesh* theMesh )
2623 const SMESHDS_Mesh* aMesh = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
2627 int nbGrp = aMesh->GetNbGroups();
2631 // iterates on groups and find necessary elements ids
2632 const std::set<SMESHDS_GroupBase*>& aGroups = aMesh->GetGroups();
2633 set<SMESHDS_GroupBase*>::const_iterator GrIt = aGroups.begin();
2634 for (; GrIt != aGroups.end(); GrIt++) {
2635 SMESHDS_GroupBase* aGrp = (*GrIt);
2638 // check type and color of group
2639 if ( !isEqual( myColor, aGrp->GetColor() ) )
2641 if ( myType != SMDSAbs_All && myType != (SMDSAbs_ElementType)aGrp->GetType() )
2644 SMDSAbs_ElementType aGrpElType = (SMDSAbs_ElementType)aGrp->GetType();
2645 if ( myType == aGrpElType || (myType == SMDSAbs_All && aGrpElType != SMDSAbs_Node) ) {
2646 // add elements IDS into control
2647 int aSize = aGrp->Extent();
2648 for (int i = 0; i < aSize; i++)
2649 myIDs.insert( aGrp->GetID(i+1) );
2654 void GroupColor::SetColorStr( const TCollection_AsciiString& theStr )
2656 TCollection_AsciiString aStr = theStr;
2657 aStr.RemoveAll( ' ' );
2658 aStr.RemoveAll( '\t' );
2659 for ( int aPos = aStr.Search( ";;" ); aPos != -1; aPos = aStr.Search( ";;" ) )
2660 aStr.Remove( aPos, 2 );
2661 Standard_Real clr[3];
2662 clr[0] = clr[1] = clr[2] = 0.;
2663 for ( int i = 0; i < 3; i++ ) {
2664 TCollection_AsciiString tmpStr = aStr.Token( ";", i+1 );
2665 if ( !tmpStr.IsEmpty() && tmpStr.IsRealValue() )
2666 clr[i] = tmpStr.RealValue();
2668 myColor = Quantity_Color( clr[0], clr[1], clr[2], Quantity_TOC_RGB );
2671 //=======================================================================
2672 // name : GetRangeStr
2673 // Purpose : Get range as a string.
2674 // Example: "1,2,3,50-60,63,67,70-"
2675 //=======================================================================
2677 void GroupColor::GetColorStr( TCollection_AsciiString& theResStr ) const
2680 theResStr += TCollection_AsciiString( myColor.Red() );
2681 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Green() );
2682 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Blue() );
2685 //================================================================================
2687 Class : ElemGeomType
2688 Description : Predicate to check element geometry type
2690 //================================================================================
2692 ElemGeomType::ElemGeomType()
2695 myType = SMDSAbs_All;
2696 myGeomType = SMDSGeom_TRIANGLE;
2699 void ElemGeomType::SetMesh( const SMDS_Mesh* theMesh )
2704 bool ElemGeomType::IsSatisfy( long theId )
2706 if (!myMesh) return false;
2707 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2710 const SMDSAbs_ElementType anElemType = anElem->GetType();
2711 if ( myType != SMDSAbs_All && anElemType != myType )
2713 bool isOk = ( anElem->GetGeomType() == myGeomType );
2717 void ElemGeomType::SetType( SMDSAbs_ElementType theType )
2722 SMDSAbs_ElementType ElemGeomType::GetType() const
2727 void ElemGeomType::SetGeomType( SMDSAbs_GeometryType theType )
2729 myGeomType = theType;
2732 SMDSAbs_GeometryType ElemGeomType::GetGeomType() const
2737 //================================================================================
2739 Class : ElemEntityType
2740 Description : Predicate to check element entity type
2742 //================================================================================
2744 ElemEntityType::ElemEntityType():
2746 myType( SMDSAbs_All ),
2747 myEntityType( SMDSEntity_0D )
2751 void ElemEntityType::SetMesh( const SMDS_Mesh* theMesh )
2756 bool ElemEntityType::IsSatisfy( long theId )
2758 if ( !myMesh ) return false;
2759 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2761 myEntityType == anElem->GetEntityType() &&
2762 ( myType == SMDSAbs_Edge || myType == SMDSAbs_Face || myType == SMDSAbs_Volume ));
2765 void ElemEntityType::SetType( SMDSAbs_ElementType theType )
2770 SMDSAbs_ElementType ElemEntityType::GetType() const
2775 void ElemEntityType::SetElemEntityType( SMDSAbs_EntityType theEntityType )
2777 myEntityType = theEntityType;
2780 SMDSAbs_EntityType ElemEntityType::GetElemEntityType() const
2782 return myEntityType;
2785 //================================================================================
2787 * \brief Class CoplanarFaces
2789 //================================================================================
2791 CoplanarFaces::CoplanarFaces()
2792 : myFaceID(0), myToler(0)
2795 void CoplanarFaces::SetMesh( const SMDS_Mesh* theMesh )
2797 myMeshModifTracer.SetMesh( theMesh );
2798 if ( myMeshModifTracer.IsMeshModified() )
2800 // Build a set of coplanar face ids
2802 myCoplanarIDs.clear();
2804 if ( !myMeshModifTracer.GetMesh() || !myFaceID || !myToler )
2807 const SMDS_MeshElement* face = myMeshModifTracer.GetMesh()->FindElement( myFaceID );
2808 if ( !face || face->GetType() != SMDSAbs_Face )
2812 gp_Vec myNorm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
2816 const double radianTol = myToler * M_PI / 180.;
2817 std::set< SMESH_TLink > checkedLinks;
2819 std::list< pair< const SMDS_MeshElement*, gp_Vec > > faceQueue;
2820 faceQueue.push_back( make_pair( face, myNorm ));
2821 while ( !faceQueue.empty() )
2823 face = faceQueue.front().first;
2824 myNorm = faceQueue.front().second;
2825 faceQueue.pop_front();
2827 for ( int i = 0, nbN = face->NbCornerNodes(); i < nbN; ++i )
2829 const SMDS_MeshNode* n1 = face->GetNode( i );
2830 const SMDS_MeshNode* n2 = face->GetNode(( i+1 )%nbN);
2831 if ( !checkedLinks.insert( SMESH_TLink( n1, n2 )).second )
2833 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator(SMDSAbs_Face);
2834 while ( fIt->more() )
2836 const SMDS_MeshElement* f = fIt->next();
2837 if ( f->GetNodeIndex( n2 ) > -1 )
2839 gp_Vec norm = getNormale( static_cast<const SMDS_MeshFace*>(f), &normOK );
2840 if (!normOK || myNorm.Angle( norm ) <= radianTol)
2842 myCoplanarIDs.insert( f->GetID() );
2843 faceQueue.push_back( make_pair( f, norm ));
2851 bool CoplanarFaces::IsSatisfy( long theElementId )
2853 return myCoplanarIDs.count( theElementId );
2858 *Description : Predicate for Range of Ids.
2859 * Range may be specified with two ways.
2860 * 1. Using AddToRange method
2861 * 2. With SetRangeStr method. Parameter of this method is a string
2862 * like as "1,2,3,50-60,63,67,70-"
2865 //=======================================================================
2866 // name : RangeOfIds
2867 // Purpose : Constructor
2868 //=======================================================================
2869 RangeOfIds::RangeOfIds()
2872 myType = SMDSAbs_All;
2875 //=======================================================================
2877 // Purpose : Set mesh
2878 //=======================================================================
2879 void RangeOfIds::SetMesh( const SMDS_Mesh* theMesh )
2884 //=======================================================================
2885 // name : AddToRange
2886 // Purpose : Add ID to the range
2887 //=======================================================================
2888 bool RangeOfIds::AddToRange( long theEntityId )
2890 myIds.Add( theEntityId );
2894 //=======================================================================
2895 // name : GetRangeStr
2896 // Purpose : Get range as a string.
2897 // Example: "1,2,3,50-60,63,67,70-"
2898 //=======================================================================
2899 void RangeOfIds::GetRangeStr( TCollection_AsciiString& theResStr )
2903 TColStd_SequenceOfInteger anIntSeq;
2904 TColStd_SequenceOfAsciiString aStrSeq;
2906 TColStd_MapIteratorOfMapOfInteger anIter( myIds );
2907 for ( ; anIter.More(); anIter.Next() )
2909 int anId = anIter.Key();
2910 TCollection_AsciiString aStr( anId );
2911 anIntSeq.Append( anId );
2912 aStrSeq.Append( aStr );
2915 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2917 int aMinId = myMin( i );
2918 int aMaxId = myMax( i );
2920 TCollection_AsciiString aStr;
2921 if ( aMinId != IntegerFirst() )
2926 if ( aMaxId != IntegerLast() )
2929 // find position of the string in result sequence and insert string in it
2930 if ( anIntSeq.Length() == 0 )
2932 anIntSeq.Append( aMinId );
2933 aStrSeq.Append( aStr );
2937 if ( aMinId < anIntSeq.First() )
2939 anIntSeq.Prepend( aMinId );
2940 aStrSeq.Prepend( aStr );
2942 else if ( aMinId > anIntSeq.Last() )
2944 anIntSeq.Append( aMinId );
2945 aStrSeq.Append( aStr );
2948 for ( int j = 1, k = anIntSeq.Length(); j <= k; j++ )
2949 if ( aMinId < anIntSeq( j ) )
2951 anIntSeq.InsertBefore( j, aMinId );
2952 aStrSeq.InsertBefore( j, aStr );
2958 if ( aStrSeq.Length() == 0 )
2961 theResStr = aStrSeq( 1 );
2962 for ( int j = 2, k = aStrSeq.Length(); j <= k; j++ )
2965 theResStr += aStrSeq( j );
2969 //=======================================================================
2970 // name : SetRangeStr
2971 // Purpose : Define range with string
2972 // Example of entry string: "1,2,3,50-60,63,67,70-"
2973 //=======================================================================
2974 bool RangeOfIds::SetRangeStr( const TCollection_AsciiString& theStr )
2980 TCollection_AsciiString aStr = theStr;
2981 aStr.RemoveAll( ' ' );
2982 aStr.RemoveAll( '\t' );
2984 for ( int aPos = aStr.Search( ",," ); aPos != -1; aPos = aStr.Search( ",," ) )
2985 aStr.Remove( aPos, 2 );
2987 TCollection_AsciiString tmpStr = aStr.Token( ",", 1 );
2989 while ( tmpStr != "" )
2991 tmpStr = aStr.Token( ",", i++ );
2992 int aPos = tmpStr.Search( '-' );
2996 if ( tmpStr.IsIntegerValue() )
2997 myIds.Add( tmpStr.IntegerValue() );
3003 TCollection_AsciiString aMaxStr = tmpStr.Split( aPos );
3004 TCollection_AsciiString aMinStr = tmpStr;
3006 while ( aMinStr.Search( "-" ) != -1 ) aMinStr.RemoveAll( '-' );
3007 while ( aMaxStr.Search( "-" ) != -1 ) aMaxStr.RemoveAll( '-' );
3009 if ( (!aMinStr.IsEmpty() && !aMinStr.IsIntegerValue()) ||
3010 (!aMaxStr.IsEmpty() && !aMaxStr.IsIntegerValue()) )
3013 myMin.Append( aMinStr.IsEmpty() ? IntegerFirst() : aMinStr.IntegerValue() );
3014 myMax.Append( aMaxStr.IsEmpty() ? IntegerLast() : aMaxStr.IntegerValue() );
3021 //=======================================================================
3023 // Purpose : Get type of supported entities
3024 //=======================================================================
3025 SMDSAbs_ElementType RangeOfIds::GetType() const
3030 //=======================================================================
3032 // Purpose : Set type of supported entities
3033 //=======================================================================
3034 void RangeOfIds::SetType( SMDSAbs_ElementType theType )
3039 //=======================================================================
3041 // Purpose : Verify whether entity satisfies to this rpedicate
3042 //=======================================================================
3043 bool RangeOfIds::IsSatisfy( long theId )
3048 if ( myType == SMDSAbs_Node )
3050 if ( myMesh->FindNode( theId ) == 0 )
3055 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
3056 if ( anElem == 0 || (myType != anElem->GetType() && myType != SMDSAbs_All ))
3060 if ( myIds.Contains( theId ) )
3063 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
3064 if ( theId >= myMin( i ) && theId <= myMax( i ) )
3072 Description : Base class for comparators
3074 Comparator::Comparator():
3078 Comparator::~Comparator()
3081 void Comparator::SetMesh( const SMDS_Mesh* theMesh )
3084 myFunctor->SetMesh( theMesh );
3087 void Comparator::SetMargin( double theValue )
3089 myMargin = theValue;
3092 void Comparator::SetNumFunctor( NumericalFunctorPtr theFunct )
3094 myFunctor = theFunct;
3097 SMDSAbs_ElementType Comparator::GetType() const
3099 return myFunctor ? myFunctor->GetType() : SMDSAbs_All;
3102 double Comparator::GetMargin()
3110 Description : Comparator "<"
3112 bool LessThan::IsSatisfy( long theId )
3114 return myFunctor && myFunctor->GetValue( theId ) < myMargin;
3120 Description : Comparator ">"
3122 bool MoreThan::IsSatisfy( long theId )
3124 return myFunctor && myFunctor->GetValue( theId ) > myMargin;
3130 Description : Comparator "="
3133 myToler(Precision::Confusion())
3136 bool EqualTo::IsSatisfy( long theId )
3138 return myFunctor && fabs( myFunctor->GetValue( theId ) - myMargin ) < myToler;
3141 void EqualTo::SetTolerance( double theToler )
3146 double EqualTo::GetTolerance()
3153 Description : Logical NOT predicate
3155 LogicalNOT::LogicalNOT()
3158 LogicalNOT::~LogicalNOT()
3161 bool LogicalNOT::IsSatisfy( long theId )
3163 return myPredicate && !myPredicate->IsSatisfy( theId );
3166 void LogicalNOT::SetMesh( const SMDS_Mesh* theMesh )
3169 myPredicate->SetMesh( theMesh );
3172 void LogicalNOT::SetPredicate( PredicatePtr thePred )
3174 myPredicate = thePred;
3177 SMDSAbs_ElementType LogicalNOT::GetType() const
3179 return myPredicate ? myPredicate->GetType() : SMDSAbs_All;
3184 Class : LogicalBinary
3185 Description : Base class for binary logical predicate
3187 LogicalBinary::LogicalBinary()
3190 LogicalBinary::~LogicalBinary()
3193 void LogicalBinary::SetMesh( const SMDS_Mesh* theMesh )
3196 myPredicate1->SetMesh( theMesh );
3199 myPredicate2->SetMesh( theMesh );
3202 void LogicalBinary::SetPredicate1( PredicatePtr thePredicate )
3204 myPredicate1 = thePredicate;
3207 void LogicalBinary::SetPredicate2( PredicatePtr thePredicate )
3209 myPredicate2 = thePredicate;
3212 SMDSAbs_ElementType LogicalBinary::GetType() const
3214 if ( !myPredicate1 || !myPredicate2 )
3217 SMDSAbs_ElementType aType1 = myPredicate1->GetType();
3218 SMDSAbs_ElementType aType2 = myPredicate2->GetType();
3220 return aType1 == aType2 ? aType1 : SMDSAbs_All;
3226 Description : Logical AND
3228 bool LogicalAND::IsSatisfy( long theId )
3233 myPredicate1->IsSatisfy( theId ) &&
3234 myPredicate2->IsSatisfy( theId );
3240 Description : Logical OR
3242 bool LogicalOR::IsSatisfy( long theId )
3247 (myPredicate1->IsSatisfy( theId ) ||
3248 myPredicate2->IsSatisfy( theId ));
3262 void Filter::SetPredicate( PredicatePtr thePredicate )
3264 myPredicate = thePredicate;
3267 void Filter::GetElementsId( const SMDS_Mesh* theMesh,
3268 PredicatePtr thePredicate,
3269 TIdSequence& theSequence )
3271 theSequence.clear();
3273 if ( !theMesh || !thePredicate )
3276 thePredicate->SetMesh( theMesh );
3278 SMDS_ElemIteratorPtr elemIt = theMesh->elementsIterator( thePredicate->GetType() );
3280 while ( elemIt->more() ) {
3281 const SMDS_MeshElement* anElem = elemIt->next();
3282 long anId = anElem->GetID();
3283 if ( thePredicate->IsSatisfy( anId ) )
3284 theSequence.push_back( anId );
3289 void Filter::GetElementsId( const SMDS_Mesh* theMesh,
3290 Filter::TIdSequence& theSequence )
3292 GetElementsId(theMesh,myPredicate,theSequence);
3299 typedef std::set<SMDS_MeshFace*> TMapOfFacePtr;
3305 ManifoldPart::Link::Link( SMDS_MeshNode* theNode1,
3306 SMDS_MeshNode* theNode2 )
3312 ManifoldPart::Link::~Link()
3318 bool ManifoldPart::Link::IsEqual( const ManifoldPart::Link& theLink ) const
3320 if ( myNode1 == theLink.myNode1 &&
3321 myNode2 == theLink.myNode2 )
3323 else if ( myNode1 == theLink.myNode2 &&
3324 myNode2 == theLink.myNode1 )
3330 bool ManifoldPart::Link::operator<( const ManifoldPart::Link& x ) const
3332 if(myNode1 < x.myNode1) return true;
3333 if(myNode1 == x.myNode1)
3334 if(myNode2 < x.myNode2) return true;
3338 bool ManifoldPart::IsEqual( const ManifoldPart::Link& theLink1,
3339 const ManifoldPart::Link& theLink2 )
3341 return theLink1.IsEqual( theLink2 );
3344 ManifoldPart::ManifoldPart()
3347 myAngToler = Precision::Angular();
3348 myIsOnlyManifold = true;
3351 ManifoldPart::~ManifoldPart()
3356 void ManifoldPart::SetMesh( const SMDS_Mesh* theMesh )
3362 SMDSAbs_ElementType ManifoldPart::GetType() const
3363 { return SMDSAbs_Face; }
3365 bool ManifoldPart::IsSatisfy( long theElementId )
3367 return myMapIds.Contains( theElementId );
3370 void ManifoldPart::SetAngleTolerance( const double theAngToler )
3371 { myAngToler = theAngToler; }
3373 double ManifoldPart::GetAngleTolerance() const
3374 { return myAngToler; }
3376 void ManifoldPart::SetIsOnlyManifold( const bool theIsOnly )
3377 { myIsOnlyManifold = theIsOnly; }
3379 void ManifoldPart::SetStartElem( const long theStartId )
3380 { myStartElemId = theStartId; }
3382 bool ManifoldPart::process()
3385 myMapBadGeomIds.Clear();
3387 myAllFacePtr.clear();
3388 myAllFacePtrIntDMap.clear();
3392 // collect all faces into own map
3393 SMDS_FaceIteratorPtr anFaceItr = myMesh->facesIterator();
3394 for (; anFaceItr->more(); )
3396 SMDS_MeshFace* aFacePtr = (SMDS_MeshFace*)anFaceItr->next();
3397 myAllFacePtr.push_back( aFacePtr );
3398 myAllFacePtrIntDMap[aFacePtr] = myAllFacePtr.size()-1;
3401 SMDS_MeshFace* aStartFace = (SMDS_MeshFace*)myMesh->FindElement( myStartElemId );
3405 // the map of non manifold links and bad geometry
3406 TMapOfLink aMapOfNonManifold;
3407 TColStd_MapOfInteger aMapOfTreated;
3409 // begin cycle on faces from start index and run on vector till the end
3410 // and from begin to start index to cover whole vector
3411 const int aStartIndx = myAllFacePtrIntDMap[aStartFace];
3412 bool isStartTreat = false;
3413 for ( int fi = aStartIndx; !isStartTreat || fi != aStartIndx ; fi++ )
3415 if ( fi == aStartIndx )
3416 isStartTreat = true;
3417 // as result next time when fi will be equal to aStartIndx
3419 SMDS_MeshFace* aFacePtr = myAllFacePtr[ fi ];
3420 if ( aMapOfTreated.Contains( aFacePtr->GetID() ) )
3423 aMapOfTreated.Add( aFacePtr->GetID() );
3424 TColStd_MapOfInteger aResFaces;
3425 if ( !findConnected( myAllFacePtrIntDMap, aFacePtr,
3426 aMapOfNonManifold, aResFaces ) )
3428 TColStd_MapIteratorOfMapOfInteger anItr( aResFaces );
3429 for ( ; anItr.More(); anItr.Next() )
3431 int aFaceId = anItr.Key();
3432 aMapOfTreated.Add( aFaceId );
3433 myMapIds.Add( aFaceId );
3436 if ( fi == ( myAllFacePtr.size() - 1 ) )
3438 } // end run on vector of faces
3439 return !myMapIds.IsEmpty();
3442 static void getLinks( const SMDS_MeshFace* theFace,
3443 ManifoldPart::TVectorOfLink& theLinks )
3445 int aNbNode = theFace->NbNodes();
3446 SMDS_ElemIteratorPtr aNodeItr = theFace->nodesIterator();
3448 SMDS_MeshNode* aNode = 0;
3449 for ( ; aNodeItr->more() && i <= aNbNode; )
3452 SMDS_MeshNode* aN1 = (SMDS_MeshNode*)aNodeItr->next();
3456 SMDS_MeshNode* aN2 = ( i >= aNbNode ) ? aNode : (SMDS_MeshNode*)aNodeItr->next();
3458 ManifoldPart::Link aLink( aN1, aN2 );
3459 theLinks.push_back( aLink );
3463 bool ManifoldPart::findConnected
3464 ( const ManifoldPart::TDataMapFacePtrInt& theAllFacePtrInt,
3465 SMDS_MeshFace* theStartFace,
3466 ManifoldPart::TMapOfLink& theNonManifold,
3467 TColStd_MapOfInteger& theResFaces )
3469 theResFaces.Clear();
3470 if ( !theAllFacePtrInt.size() )
3473 if ( getNormale( theStartFace ).SquareModulus() <= gp::Resolution() )
3475 myMapBadGeomIds.Add( theStartFace->GetID() );
3479 ManifoldPart::TMapOfLink aMapOfBoundary, aMapToSkip;
3480 ManifoldPart::TVectorOfLink aSeqOfBoundary;
3481 theResFaces.Add( theStartFace->GetID() );
3482 ManifoldPart::TDataMapOfLinkFacePtr aDMapLinkFace;
3484 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3485 aDMapLinkFace, theNonManifold, theStartFace );
3487 bool isDone = false;
3488 while ( !isDone && aMapOfBoundary.size() != 0 )
3490 bool isToReset = false;
3491 ManifoldPart::TVectorOfLink::iterator pLink = aSeqOfBoundary.begin();
3492 for ( ; !isToReset && pLink != aSeqOfBoundary.end(); ++pLink )
3494 ManifoldPart::Link aLink = *pLink;
3495 if ( aMapToSkip.find( aLink ) != aMapToSkip.end() )
3497 // each link could be treated only once
3498 aMapToSkip.insert( aLink );
3500 ManifoldPart::TVectorOfFacePtr aFaces;
3502 if ( myIsOnlyManifold &&
3503 (theNonManifold.find( aLink ) != theNonManifold.end()) )
3507 getFacesByLink( aLink, aFaces );
3508 // filter the element to keep only indicated elements
3509 ManifoldPart::TVectorOfFacePtr aFiltered;
3510 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3511 for ( ; pFace != aFaces.end(); ++pFace )
3513 SMDS_MeshFace* aFace = *pFace;
3514 if ( myAllFacePtrIntDMap.find( aFace ) != myAllFacePtrIntDMap.end() )
3515 aFiltered.push_back( aFace );
3518 if ( aFaces.size() < 2 ) // no neihgbour faces
3520 else if ( myIsOnlyManifold && aFaces.size() > 2 ) // non manifold case
3522 theNonManifold.insert( aLink );
3527 // compare normal with normals of neighbor element
3528 SMDS_MeshFace* aPrevFace = aDMapLinkFace[ aLink ];
3529 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3530 for ( ; pFace != aFaces.end(); ++pFace )
3532 SMDS_MeshFace* aNextFace = *pFace;
3533 if ( aPrevFace == aNextFace )
3535 int anNextFaceID = aNextFace->GetID();
3536 if ( myIsOnlyManifold && theResFaces.Contains( anNextFaceID ) )
3537 // should not be with non manifold restriction. probably bad topology
3539 // check if face was treated and skipped
3540 if ( myMapBadGeomIds.Contains( anNextFaceID ) ||
3541 !isInPlane( aPrevFace, aNextFace ) )
3543 // add new element to connected and extend the boundaries.
3544 theResFaces.Add( anNextFaceID );
3545 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3546 aDMapLinkFace, theNonManifold, aNextFace );
3550 isDone = !isToReset;
3553 return !theResFaces.IsEmpty();
3556 bool ManifoldPart::isInPlane( const SMDS_MeshFace* theFace1,
3557 const SMDS_MeshFace* theFace2 )
3559 gp_Dir aNorm1 = gp_Dir( getNormale( theFace1 ) );
3560 gp_XYZ aNorm2XYZ = getNormale( theFace2 );
3561 if ( aNorm2XYZ.SquareModulus() <= gp::Resolution() )
3563 myMapBadGeomIds.Add( theFace2->GetID() );
3566 if ( aNorm1.IsParallel( gp_Dir( aNorm2XYZ ), myAngToler ) )
3572 void ManifoldPart::expandBoundary
3573 ( ManifoldPart::TMapOfLink& theMapOfBoundary,
3574 ManifoldPart::TVectorOfLink& theSeqOfBoundary,
3575 ManifoldPart::TDataMapOfLinkFacePtr& theDMapLinkFacePtr,
3576 ManifoldPart::TMapOfLink& theNonManifold,
3577 SMDS_MeshFace* theNextFace ) const
3579 ManifoldPart::TVectorOfLink aLinks;
3580 getLinks( theNextFace, aLinks );
3581 int aNbLink = (int)aLinks.size();
3582 for ( int i = 0; i < aNbLink; i++ )
3584 ManifoldPart::Link aLink = aLinks[ i ];
3585 if ( myIsOnlyManifold && (theNonManifold.find( aLink ) != theNonManifold.end()) )
3587 if ( theMapOfBoundary.find( aLink ) != theMapOfBoundary.end() )
3589 if ( myIsOnlyManifold )
3591 // remove from boundary
3592 theMapOfBoundary.erase( aLink );
3593 ManifoldPart::TVectorOfLink::iterator pLink = theSeqOfBoundary.begin();
3594 for ( ; pLink != theSeqOfBoundary.end(); ++pLink )
3596 ManifoldPart::Link aBoundLink = *pLink;
3597 if ( aBoundLink.IsEqual( aLink ) )
3599 theSeqOfBoundary.erase( pLink );
3607 theMapOfBoundary.insert( aLink );
3608 theSeqOfBoundary.push_back( aLink );
3609 theDMapLinkFacePtr[ aLink ] = theNextFace;
3614 void ManifoldPart::getFacesByLink( const ManifoldPart::Link& theLink,
3615 ManifoldPart::TVectorOfFacePtr& theFaces ) const
3617 std::set<SMDS_MeshCell *> aSetOfFaces;
3618 // take all faces that shared first node
3619 SMDS_ElemIteratorPtr anItr = theLink.myNode1->facesIterator();
3620 for ( ; anItr->more(); )
3622 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3625 aSetOfFaces.insert( aFace );
3627 // take all faces that shared second node
3628 anItr = theLink.myNode2->facesIterator();
3629 // find the common part of two sets
3630 for ( ; anItr->more(); )
3632 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3633 if ( aSetOfFaces.count( aFace ) )
3634 theFaces.push_back( aFace );
3643 ElementsOnSurface::ElementsOnSurface()
3646 myType = SMDSAbs_All;
3648 myToler = Precision::Confusion();
3649 myUseBoundaries = false;
3652 ElementsOnSurface::~ElementsOnSurface()
3656 void ElementsOnSurface::SetMesh( const SMDS_Mesh* theMesh )
3658 myMeshModifTracer.SetMesh( theMesh );
3659 if ( myMeshModifTracer.IsMeshModified())
3663 bool ElementsOnSurface::IsSatisfy( long theElementId )
3665 return myIds.Contains( theElementId );
3668 SMDSAbs_ElementType ElementsOnSurface::GetType() const
3671 void ElementsOnSurface::SetTolerance( const double theToler )
3673 if ( myToler != theToler )
3678 double ElementsOnSurface::GetTolerance() const
3681 void ElementsOnSurface::SetUseBoundaries( bool theUse )
3683 if ( myUseBoundaries != theUse ) {
3684 myUseBoundaries = theUse;
3685 SetSurface( mySurf, myType );
3689 void ElementsOnSurface::SetSurface( const TopoDS_Shape& theShape,
3690 const SMDSAbs_ElementType theType )
3695 if ( theShape.IsNull() || theShape.ShapeType() != TopAbs_FACE )
3697 mySurf = TopoDS::Face( theShape );
3698 BRepAdaptor_Surface SA( mySurf, myUseBoundaries );
3700 u1 = SA.FirstUParameter(),
3701 u2 = SA.LastUParameter(),
3702 v1 = SA.FirstVParameter(),
3703 v2 = SA.LastVParameter();
3704 Handle(Geom_Surface) surf = BRep_Tool::Surface( mySurf );
3705 myProjector.Init( surf, u1,u2, v1,v2 );
3709 void ElementsOnSurface::process()
3712 if ( mySurf.IsNull() )
3715 if ( !myMeshModifTracer.GetMesh() )
3718 myIds.ReSize( myMeshModifTracer.GetMesh()->GetMeshInfo().NbElements( myType ));
3720 SMDS_ElemIteratorPtr anIter = myMeshModifTracer.GetMesh()->elementsIterator( myType );
3721 for(; anIter->more(); )
3722 process( anIter->next() );
3725 void ElementsOnSurface::process( const SMDS_MeshElement* theElemPtr )
3727 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3728 bool isSatisfy = true;
3729 for ( ; aNodeItr->more(); )
3731 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3732 if ( !isOnSurface( aNode ) )
3739 myIds.Add( theElemPtr->GetID() );
3742 bool ElementsOnSurface::isOnSurface( const SMDS_MeshNode* theNode )
3744 if ( mySurf.IsNull() )
3747 gp_Pnt aPnt( theNode->X(), theNode->Y(), theNode->Z() );
3748 // double aToler2 = myToler * myToler;
3749 // if ( mySurf->IsKind(STANDARD_TYPE(Geom_Plane)))
3751 // gp_Pln aPln = Handle(Geom_Plane)::DownCast(mySurf)->Pln();
3752 // if ( aPln.SquareDistance( aPnt ) > aToler2 )
3755 // else if ( mySurf->IsKind(STANDARD_TYPE(Geom_CylindricalSurface)))
3757 // gp_Cylinder aCyl = Handle(Geom_CylindricalSurface)::DownCast(mySurf)->Cylinder();
3758 // double aRad = aCyl.Radius();
3759 // gp_Ax3 anAxis = aCyl.Position();
3760 // gp_XYZ aLoc = aCyl.Location().XYZ();
3761 // double aXDist = anAxis.XDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3762 // double aYDist = anAxis.YDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3763 // if ( fabs(aXDist*aXDist + aYDist*aYDist - aRad*aRad) > aToler2 )
3768 myProjector.Perform( aPnt );
3769 bool isOn = ( myProjector.IsDone() && myProjector.LowerDistance() <= myToler );
3779 ElementsOnShape::ElementsOnShape()
3781 myType(SMDSAbs_All),
3782 myToler(Precision::Confusion()),
3783 myAllNodesFlag(false)
3787 ElementsOnShape::~ElementsOnShape()
3792 SMDSAbs_ElementType ElementsOnShape::GetType() const
3797 void ElementsOnShape::SetTolerance (const double theToler)
3799 if (myToler != theToler) {
3801 SetShape(myShape, myType);
3805 double ElementsOnShape::GetTolerance() const
3810 void ElementsOnShape::SetAllNodes (bool theAllNodes)
3812 myAllNodesFlag = theAllNodes;
3815 void ElementsOnShape::SetMesh (const SMDS_Mesh* theMesh)
3820 void ElementsOnShape::SetShape (const TopoDS_Shape& theShape,
3821 const SMDSAbs_ElementType theType)
3825 if ( myShape.IsNull() ) return;
3827 TopTools_IndexedMapOfShape shapesMap;
3828 TopAbs_ShapeEnum shapeTypes[4] = { TopAbs_SOLID, TopAbs_FACE, TopAbs_EDGE, TopAbs_VERTEX };
3829 TopExp_Explorer sub;
3830 for ( int i = 0; i < 4; ++i )
3832 if ( shapesMap.IsEmpty() )
3833 for ( sub.Init( myShape, shapeTypes[i] ); sub.More(); sub.Next() )
3834 shapesMap.Add( sub.Current() );
3836 for ( sub.Init( myShape, shapeTypes[i], shapeTypes[i-1] ); sub.More(); sub.Next() )
3837 shapesMap.Add( sub.Current() );
3841 myClassifiers.resize( shapesMap.Extent() );
3842 for ( int i = 0; i < shapesMap.Extent(); ++i )
3843 myClassifiers[ i ] = new TClassifier( shapesMap( i+1 ), myToler );
3846 void ElementsOnShape::clearClassifiers()
3848 for ( size_t i = 0; i < myClassifiers.size(); ++i )
3849 delete myClassifiers[ i ];
3850 myClassifiers.clear();
3853 bool ElementsOnShape::IsSatisfy (long elemId)
3855 const SMDS_MeshElement* elem =
3856 ( myType == SMDSAbs_Node ? myMesh->FindNode( elemId ) : myMesh->FindElement( elemId ));
3857 if ( !elem || myClassifiers.empty() )
3860 for ( size_t i = 0; i < myClassifiers.size(); ++i )
3862 SMDS_ElemIteratorPtr aNodeItr = elem->nodesIterator();
3863 bool isSatisfy = myAllNodesFlag;
3865 gp_XYZ centerXYZ (0, 0, 0);
3867 while (aNodeItr->more() && (isSatisfy == myAllNodesFlag))
3869 SMESH_TNodeXYZ aPnt ( aNodeItr->next() );
3871 isSatisfy = ! myClassifiers[i]->IsOut( aPnt );
3874 // Check the center point for volumes MantisBug 0020168
3877 myClassifiers[i]->ShapeType() == TopAbs_SOLID)
3879 centerXYZ /= elem->NbNodes();
3880 isSatisfy = ! myClassifiers[i]->IsOut( centerXYZ );
3889 TopAbs_ShapeEnum ElementsOnShape::TClassifier::ShapeType() const
3891 return myShape.ShapeType();
3894 bool ElementsOnShape::TClassifier::IsOut(const gp_Pnt& p)
3896 return (this->*myIsOutFun)( p );
3899 void ElementsOnShape::TClassifier::Init (const TopoDS_Shape& theShape, double theTol)
3903 switch ( myShape.ShapeType() )
3905 case TopAbs_SOLID: {
3906 mySolidClfr.Load(theShape);
3907 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfSolid;
3911 Standard_Real u1,u2,v1,v2;
3912 Handle(Geom_Surface) surf = BRep_Tool::Surface( TopoDS::Face( theShape ));
3913 surf->Bounds( u1,u2,v1,v2 );
3914 myProjFace.Init(surf, u1,u2, v1,v2, myTol );
3915 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfFace;
3919 Standard_Real u1, u2;
3920 Handle(Geom_Curve) curve = BRep_Tool::Curve( TopoDS::Edge(theShape), u1, u2);
3921 myProjEdge.Init(curve, u1, u2);
3922 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfEdge;
3925 case TopAbs_VERTEX:{
3926 myVertexXYZ = BRep_Tool::Pnt( TopoDS::Vertex( theShape ) );
3927 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfVertex;
3931 throw SALOME_Exception("Programmer error in usage of ElementsOnShape::TClassifier");
3935 bool ElementsOnShape::TClassifier::isOutOfSolid (const gp_Pnt& p)
3937 mySolidClfr.Perform( p, myTol );
3938 return ( mySolidClfr.State() != TopAbs_IN && mySolidClfr.State() != TopAbs_ON );
3941 bool ElementsOnShape::TClassifier::isOutOfFace (const gp_Pnt& p)
3943 myProjFace.Perform( p );
3944 if ( myProjFace.IsDone() && myProjFace.LowerDistance() <= myTol )
3946 // check relatively to the face
3947 Quantity_Parameter u, v;
3948 myProjFace.LowerDistanceParameters(u, v);
3949 gp_Pnt2d aProjPnt (u, v);
3950 BRepClass_FaceClassifier aClsf ( TopoDS::Face( myShape ), aProjPnt, myTol );
3951 if ( aClsf.State() == TopAbs_IN || aClsf.State() == TopAbs_ON )
3957 bool ElementsOnShape::TClassifier::isOutOfEdge (const gp_Pnt& p)
3959 myProjEdge.Perform( p );
3960 return ! ( myProjEdge.NbPoints() > 0 && myProjEdge.LowerDistance() <= myTol );
3963 bool ElementsOnShape::TClassifier::isOutOfVertex(const gp_Pnt& p)
3965 return ( myVertexXYZ.Distance( p ) > myTol );
3969 TSequenceOfXYZ::TSequenceOfXYZ()
3972 TSequenceOfXYZ::TSequenceOfXYZ(size_type n) : myArray(n)
3975 TSequenceOfXYZ::TSequenceOfXYZ(size_type n, const gp_XYZ& t) : myArray(n,t)
3978 TSequenceOfXYZ::TSequenceOfXYZ(const TSequenceOfXYZ& theSequenceOfXYZ) : myArray(theSequenceOfXYZ.myArray)
3981 template <class InputIterator>
3982 TSequenceOfXYZ::TSequenceOfXYZ(InputIterator theBegin, InputIterator theEnd): myArray(theBegin,theEnd)
3985 TSequenceOfXYZ::~TSequenceOfXYZ()
3988 TSequenceOfXYZ& TSequenceOfXYZ::operator=(const TSequenceOfXYZ& theSequenceOfXYZ)
3990 myArray = theSequenceOfXYZ.myArray;
3994 gp_XYZ& TSequenceOfXYZ::operator()(size_type n)
3996 return myArray[n-1];
3999 const gp_XYZ& TSequenceOfXYZ::operator()(size_type n) const
4001 return myArray[n-1];
4004 void TSequenceOfXYZ::clear()
4009 void TSequenceOfXYZ::reserve(size_type n)
4014 void TSequenceOfXYZ::push_back(const gp_XYZ& v)
4016 myArray.push_back(v);
4019 TSequenceOfXYZ::size_type TSequenceOfXYZ::size() const
4021 return myArray.size();
4024 TMeshModifTracer::TMeshModifTracer():
4025 myMeshModifTime(0), myMesh(0)
4028 void TMeshModifTracer::SetMesh( const SMDS_Mesh* theMesh )
4030 if ( theMesh != myMesh )
4031 myMeshModifTime = 0;
4034 bool TMeshModifTracer::IsMeshModified()
4036 bool modified = false;
4039 modified = ( myMeshModifTime != myMesh->GetMTime() );
4040 myMeshModifTime = myMesh->GetMTime();