1 // Copyright (C) 2007-2015 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, or (at your option) any later version.
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
23 #include "SMESH_ControlsDef.hxx"
25 #include "SMDS_BallElement.hxx"
26 #include "SMDS_Iterator.hxx"
27 #include "SMDS_Mesh.hxx"
28 #include "SMDS_MeshElement.hxx"
29 #include "SMDS_MeshNode.hxx"
30 #include "SMDS_QuadraticEdge.hxx"
31 #include "SMDS_QuadraticFaceOfNodes.hxx"
32 #include "SMDS_VolumeTool.hxx"
33 #include "SMESHDS_GroupBase.hxx"
34 #include "SMESHDS_Mesh.hxx"
35 #include "SMESH_OctreeNode.hxx"
36 #include "SMESH_MeshAlgos.hxx"
38 #include <Basics_Utils.hxx>
40 #include <BRepAdaptor_Surface.hxx>
41 #include <BRepClass_FaceClassifier.hxx>
42 #include <BRep_Tool.hxx>
43 #include <Geom_CylindricalSurface.hxx>
44 #include <Geom_Plane.hxx>
45 #include <Geom_Surface.hxx>
46 #include <Precision.hxx>
47 #include <TColStd_MapIteratorOfMapOfInteger.hxx>
48 #include <TColStd_MapOfInteger.hxx>
49 #include <TColStd_SequenceOfAsciiString.hxx>
50 #include <TColgp_Array1OfXYZ.hxx>
54 #include <TopoDS_Edge.hxx>
55 #include <TopoDS_Face.hxx>
56 #include <TopoDS_Iterator.hxx>
57 #include <TopoDS_Shape.hxx>
58 #include <TopoDS_Vertex.hxx>
60 #include <gp_Cylinder.hxx>
67 #include <vtkMeshQuality.h>
78 const double theEps = 1e-100;
79 const double theInf = 1e+100;
81 inline gp_XYZ gpXYZ(const SMDS_MeshNode* aNode )
83 return gp_XYZ(aNode->X(), aNode->Y(), aNode->Z() );
86 inline double getAngle( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
88 gp_Vec v1( P1 - P2 ), v2( P3 - P2 );
90 return v1.Magnitude() < gp::Resolution() ||
91 v2.Magnitude() < gp::Resolution() ? 0 : v1.Angle( v2 );
94 inline double getArea( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
96 gp_Vec aVec1( P2 - P1 );
97 gp_Vec aVec2( P3 - P1 );
98 return ( aVec1 ^ aVec2 ).Magnitude() * 0.5;
101 inline double getArea( const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& P3 )
103 return getArea( P1.XYZ(), P2.XYZ(), P3.XYZ() );
108 inline double getDistance( const gp_XYZ& P1, const gp_XYZ& P2 )
110 double aDist = gp_Pnt( P1 ).Distance( gp_Pnt( P2 ) );
114 int getNbMultiConnection( const SMDS_Mesh* theMesh, const int theId )
119 const SMDS_MeshElement* anEdge = theMesh->FindElement( theId );
120 if ( anEdge == 0 || anEdge->GetType() != SMDSAbs_Edge/* || anEdge->NbNodes() != 2 */)
123 // for each pair of nodes in anEdge (there are 2 pairs in a quadratic edge)
124 // count elements containing both nodes of the pair.
125 // Note that there may be such cases for a quadratic edge (a horizontal line):
130 // +-----+------+ +-----+------+
133 // result sould be 2 in both cases
135 int aResult0 = 0, aResult1 = 0;
136 // last node, it is a medium one in a quadratic edge
137 const SMDS_MeshNode* aLastNode = anEdge->GetNode( anEdge->NbNodes() - 1 );
138 const SMDS_MeshNode* aNode0 = anEdge->GetNode( 0 );
139 const SMDS_MeshNode* aNode1 = anEdge->GetNode( 1 );
140 if ( aNode1 == aLastNode ) aNode1 = 0;
142 SMDS_ElemIteratorPtr anElemIter = aLastNode->GetInverseElementIterator();
143 while( anElemIter->more() ) {
144 const SMDS_MeshElement* anElem = anElemIter->next();
145 if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
146 SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
147 while ( anIter->more() ) {
148 if ( const SMDS_MeshElement* anElemNode = anIter->next() ) {
149 if ( anElemNode == aNode0 ) {
151 if ( !aNode1 ) break; // not a quadratic edge
153 else if ( anElemNode == aNode1 )
159 int aResult = std::max ( aResult0, aResult1 );
164 gp_XYZ getNormale( const SMDS_MeshFace* theFace, bool* ok=0 )
166 int aNbNode = theFace->NbNodes();
168 gp_XYZ q1 = gpXYZ( theFace->GetNode(1)) - gpXYZ( theFace->GetNode(0));
169 gp_XYZ q2 = gpXYZ( theFace->GetNode(2)) - gpXYZ( theFace->GetNode(0));
172 gp_XYZ q3 = gpXYZ( theFace->GetNode(3)) - gpXYZ( theFace->GetNode(0));
175 double len = n.Modulus();
176 bool zeroLen = ( len <= numeric_limits<double>::min());
180 if (ok) *ok = !zeroLen;
188 using namespace SMESH::Controls;
194 //================================================================================
196 Class : NumericalFunctor
197 Description : Base class for numerical functors
199 //================================================================================
201 NumericalFunctor::NumericalFunctor():
207 void NumericalFunctor::SetMesh( const SMDS_Mesh* theMesh )
212 bool NumericalFunctor::GetPoints(const int theId,
213 TSequenceOfXYZ& theRes ) const
220 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
221 if ( !anElem || anElem->GetType() != this->GetType() )
224 return GetPoints( anElem, theRes );
227 bool NumericalFunctor::GetPoints(const SMDS_MeshElement* anElem,
228 TSequenceOfXYZ& theRes )
235 theRes.reserve( anElem->NbNodes() );
236 theRes.setElement( anElem );
238 // Get nodes of the element
239 SMDS_ElemIteratorPtr anIter;
241 if ( anElem->IsQuadratic() ) {
242 switch ( anElem->GetType() ) {
244 anIter = dynamic_cast<const SMDS_VtkEdge*>
245 (anElem)->interlacedNodesElemIterator();
248 anIter = dynamic_cast<const SMDS_VtkFace*>
249 (anElem)->interlacedNodesElemIterator();
252 anIter = anElem->nodesIterator();
256 anIter = anElem->nodesIterator();
261 while( anIter->more() ) {
262 if ( const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>( anIter->next() ))
264 aNode->GetXYZ( xyz );
265 theRes.push_back( gp_XYZ( xyz[0], xyz[1], xyz[2] ));
273 long NumericalFunctor::GetPrecision() const
278 void NumericalFunctor::SetPrecision( const long thePrecision )
280 myPrecision = thePrecision;
281 myPrecisionValue = pow( 10., (double)( myPrecision ) );
284 double NumericalFunctor::GetValue( long theId )
288 myCurrElement = myMesh->FindElement( theId );
291 if ( GetPoints( theId, P )) // elem type is checked here
292 aVal = Round( GetValue( P ));
297 double NumericalFunctor::Round( const double & aVal )
299 return ( myPrecision >= 0 ) ? floor( aVal * myPrecisionValue + 0.5 ) / myPrecisionValue : aVal;
302 //================================================================================
304 * \brief Return histogram of functor values
305 * \param nbIntervals - number of intervals
306 * \param nbEvents - number of mesh elements having values within i-th interval
307 * \param funValues - boundaries of intervals
308 * \param elements - elements to check vulue of; empty list means "of all"
309 * \param minmax - boundaries of diapason of values to divide into intervals
311 //================================================================================
313 void NumericalFunctor::GetHistogram(int nbIntervals,
314 std::vector<int>& nbEvents,
315 std::vector<double>& funValues,
316 const vector<int>& elements,
317 const double* minmax,
318 const bool isLogarithmic)
320 if ( nbIntervals < 1 ||
322 !myMesh->GetMeshInfo().NbElements( GetType() ))
324 nbEvents.resize( nbIntervals, 0 );
325 funValues.resize( nbIntervals+1 );
327 // get all values sorted
328 std::multiset< double > values;
329 if ( elements.empty() )
331 SMDS_ElemIteratorPtr elemIt = myMesh->elementsIterator( GetType() );
332 while ( elemIt->more() )
333 values.insert( GetValue( elemIt->next()->GetID() ));
337 vector<int>::const_iterator id = elements.begin();
338 for ( ; id != elements.end(); ++id )
339 values.insert( GetValue( *id ));
344 funValues[0] = minmax[0];
345 funValues[nbIntervals] = minmax[1];
349 funValues[0] = *values.begin();
350 funValues[nbIntervals] = *values.rbegin();
352 // case nbIntervals == 1
353 if ( nbIntervals == 1 )
355 nbEvents[0] = values.size();
359 if (funValues.front() == funValues.back())
361 nbEvents.resize( 1 );
362 nbEvents[0] = values.size();
363 funValues[1] = funValues.back();
364 funValues.resize( 2 );
367 std::multiset< double >::iterator min = values.begin(), max;
368 for ( int i = 0; i < nbIntervals; ++i )
370 // find end value of i-th interval
371 double r = (i+1) / double(nbIntervals);
372 if (isLogarithmic && funValues.front() > 1e-07 && funValues.back() > 1e-07) {
373 double logmin = log10(funValues.front());
374 double lval = logmin + r * (log10(funValues.back()) - logmin);
375 funValues[i+1] = pow(10.0, lval);
378 funValues[i+1] = funValues.front() * (1-r) + funValues.back() * r;
381 // count values in the i-th interval if there are any
382 if ( min != values.end() && *min <= funValues[i+1] )
384 // find the first value out of the interval
385 max = values.upper_bound( funValues[i+1] ); // max is greater than funValues[i+1], or end()
386 nbEvents[i] = std::distance( min, max );
390 // add values larger than minmax[1]
391 nbEvents.back() += std::distance( min, values.end() );
394 //=======================================================================
397 Description : Functor calculating volume of a 3D element
399 //================================================================================
401 double Volume::GetValue( long theElementId )
403 if ( theElementId && myMesh ) {
404 SMDS_VolumeTool aVolumeTool;
405 if ( aVolumeTool.Set( myMesh->FindElement( theElementId )))
406 return aVolumeTool.GetSize();
411 double Volume::GetBadRate( double Value, int /*nbNodes*/ ) const
416 SMDSAbs_ElementType Volume::GetType() const
418 return SMDSAbs_Volume;
421 //=======================================================================
423 Class : MaxElementLength2D
424 Description : Functor calculating maximum length of 2D element
426 //================================================================================
428 double MaxElementLength2D::GetValue( const TSequenceOfXYZ& P )
434 if( len == 3 ) { // triangles
435 double L1 = getDistance(P( 1 ),P( 2 ));
436 double L2 = getDistance(P( 2 ),P( 3 ));
437 double L3 = getDistance(P( 3 ),P( 1 ));
438 aVal = Max(L1,Max(L2,L3));
440 else if( len == 4 ) { // quadrangles
441 double L1 = getDistance(P( 1 ),P( 2 ));
442 double L2 = getDistance(P( 2 ),P( 3 ));
443 double L3 = getDistance(P( 3 ),P( 4 ));
444 double L4 = getDistance(P( 4 ),P( 1 ));
445 double D1 = getDistance(P( 1 ),P( 3 ));
446 double D2 = getDistance(P( 2 ),P( 4 ));
447 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
449 else if( len == 6 ) { // quadratic triangles
450 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
451 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
452 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
453 aVal = Max(L1,Max(L2,L3));
455 else if( len == 8 || len == 9 ) { // quadratic quadrangles
456 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
457 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
458 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
459 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
460 double D1 = getDistance(P( 1 ),P( 5 ));
461 double D2 = getDistance(P( 3 ),P( 7 ));
462 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
464 // Diagonals are undefined for concave polygons
465 // else if ( P.getElementEntity() == SMDSEntity_Quad_Polygon && P.size() > 2 ) // quad polygon
468 // aVal = getDistance( P( 1 ), P( P.size() )) + getDistance( P( P.size() ), P( P.size()-1 ));
469 // for ( size_t i = 1; i < P.size()-1; i += 2 )
471 // double L = getDistance( P( i ), P( i+1 )) + getDistance( P( i+1 ), P( i+2 ));
472 // aVal = Max( aVal, L );
475 // for ( int i = P.size()-5; i > 0; i -= 2 )
476 // for ( int j = i + 4; j < P.size() + i - 2; i += 2 )
478 // double D = getDistance( P( i ), P( j ));
479 // aVal = Max( aVal, D );
486 if( myPrecision >= 0 )
488 double prec = pow( 10., (double)myPrecision );
489 aVal = floor( aVal * prec + 0.5 ) / prec;
494 double MaxElementLength2D::GetValue( long theElementId )
497 return GetPoints( theElementId, P ) ? GetValue(P) : 0.0;
500 double MaxElementLength2D::GetBadRate( double Value, int /*nbNodes*/ ) const
505 SMDSAbs_ElementType MaxElementLength2D::GetType() const
510 //=======================================================================
512 Class : MaxElementLength3D
513 Description : Functor calculating maximum length of 3D element
515 //================================================================================
517 double MaxElementLength3D::GetValue( long theElementId )
520 if( GetPoints( theElementId, P ) ) {
522 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
523 SMDSAbs_ElementType aType = aElem->GetType();
527 if( len == 4 ) { // tetras
528 double L1 = getDistance(P( 1 ),P( 2 ));
529 double L2 = getDistance(P( 2 ),P( 3 ));
530 double L3 = getDistance(P( 3 ),P( 1 ));
531 double L4 = getDistance(P( 1 ),P( 4 ));
532 double L5 = getDistance(P( 2 ),P( 4 ));
533 double L6 = getDistance(P( 3 ),P( 4 ));
534 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
537 else if( len == 5 ) { // pyramids
538 double L1 = getDistance(P( 1 ),P( 2 ));
539 double L2 = getDistance(P( 2 ),P( 3 ));
540 double L3 = getDistance(P( 3 ),P( 4 ));
541 double L4 = getDistance(P( 4 ),P( 1 ));
542 double L5 = getDistance(P( 1 ),P( 5 ));
543 double L6 = getDistance(P( 2 ),P( 5 ));
544 double L7 = getDistance(P( 3 ),P( 5 ));
545 double L8 = getDistance(P( 4 ),P( 5 ));
546 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
547 aVal = Max(aVal,Max(L7,L8));
550 else if( len == 6 ) { // pentas
551 double L1 = getDistance(P( 1 ),P( 2 ));
552 double L2 = getDistance(P( 2 ),P( 3 ));
553 double L3 = getDistance(P( 3 ),P( 1 ));
554 double L4 = getDistance(P( 4 ),P( 5 ));
555 double L5 = getDistance(P( 5 ),P( 6 ));
556 double L6 = getDistance(P( 6 ),P( 4 ));
557 double L7 = getDistance(P( 1 ),P( 4 ));
558 double L8 = getDistance(P( 2 ),P( 5 ));
559 double L9 = getDistance(P( 3 ),P( 6 ));
560 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
561 aVal = Max(aVal,Max(Max(L7,L8),L9));
564 else if( len == 8 ) { // hexas
565 double L1 = getDistance(P( 1 ),P( 2 ));
566 double L2 = getDistance(P( 2 ),P( 3 ));
567 double L3 = getDistance(P( 3 ),P( 4 ));
568 double L4 = getDistance(P( 4 ),P( 1 ));
569 double L5 = getDistance(P( 5 ),P( 6 ));
570 double L6 = getDistance(P( 6 ),P( 7 ));
571 double L7 = getDistance(P( 7 ),P( 8 ));
572 double L8 = getDistance(P( 8 ),P( 5 ));
573 double L9 = getDistance(P( 1 ),P( 5 ));
574 double L10= getDistance(P( 2 ),P( 6 ));
575 double L11= getDistance(P( 3 ),P( 7 ));
576 double L12= getDistance(P( 4 ),P( 8 ));
577 double D1 = getDistance(P( 1 ),P( 7 ));
578 double D2 = getDistance(P( 2 ),P( 8 ));
579 double D3 = getDistance(P( 3 ),P( 5 ));
580 double D4 = getDistance(P( 4 ),P( 6 ));
581 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
582 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
583 aVal = Max(aVal,Max(L11,L12));
584 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
587 else if( len == 12 ) { // hexagonal prism
588 for ( int i1 = 1; i1 < 12; ++i1 )
589 for ( int i2 = i1+1; i1 <= 12; ++i1 )
590 aVal = Max( aVal, getDistance(P( i1 ),P( i2 )));
593 else if( len == 10 ) { // quadratic tetras
594 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
595 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
596 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
597 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
598 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
599 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
600 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
603 else if( len == 13 ) { // quadratic pyramids
604 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
605 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
606 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
607 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
608 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
609 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
610 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
611 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
612 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
613 aVal = Max(aVal,Max(L7,L8));
616 else if( len == 15 ) { // quadratic pentas
617 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
618 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
619 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
620 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
621 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
622 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
623 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
624 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
625 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
626 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
627 aVal = Max(aVal,Max(Max(L7,L8),L9));
630 else if( len == 20 || len == 27 ) { // quadratic hexas
631 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
632 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
633 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
634 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
635 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
636 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
637 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
638 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
639 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
640 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
641 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
642 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
643 double D1 = getDistance(P( 1 ),P( 7 ));
644 double D2 = getDistance(P( 2 ),P( 8 ));
645 double D3 = getDistance(P( 3 ),P( 5 ));
646 double D4 = getDistance(P( 4 ),P( 6 ));
647 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
648 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
649 aVal = Max(aVal,Max(L11,L12));
650 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
653 else if( len > 1 && aElem->IsPoly() ) { // polys
654 // get the maximum distance between all pairs of nodes
655 for( int i = 1; i <= len; i++ ) {
656 for( int j = 1; j <= len; j++ ) {
657 if( j > i ) { // optimization of the loop
658 double D = getDistance( P(i), P(j) );
659 aVal = Max( aVal, D );
666 if( myPrecision >= 0 )
668 double prec = pow( 10., (double)myPrecision );
669 aVal = floor( aVal * prec + 0.5 ) / prec;
676 double MaxElementLength3D::GetBadRate( double Value, int /*nbNodes*/ ) const
681 SMDSAbs_ElementType MaxElementLength3D::GetType() const
683 return SMDSAbs_Volume;
686 //=======================================================================
689 Description : Functor for calculation of minimum angle
691 //================================================================================
693 double MinimumAngle::GetValue( const TSequenceOfXYZ& P )
700 aMin = getAngle(P( P.size() ), P( 1 ), P( 2 ));
701 aMin = Min(aMin,getAngle(P( P.size()-1 ), P( P.size() ), P( 1 )));
703 for ( int i = 2; i < P.size(); i++ )
705 double A0 = getAngle( P( i-1 ), P( i ), P( i+1 ) );
709 return aMin * 180.0 / M_PI;
712 double MinimumAngle::GetBadRate( double Value, int nbNodes ) const
714 //const double aBestAngle = PI / nbNodes;
715 const double aBestAngle = 180.0 - ( 360.0 / double(nbNodes) );
716 return ( fabs( aBestAngle - Value ));
719 SMDSAbs_ElementType MinimumAngle::GetType() const
725 //================================================================================
728 Description : Functor for calculating aspect ratio
730 //================================================================================
732 double AspectRatio::GetValue( long theId )
735 myCurrElement = myMesh->FindElement( theId );
736 if ( myCurrElement && myCurrElement->GetVtkType() == VTK_QUAD )
739 vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myCurrElement->getMeshId()]->getGrid();
740 if ( vtkCell* avtkCell = grid->GetCell( myCurrElement->getVtkId() ))
741 aVal = Round( vtkMeshQuality::QuadAspectRatio( avtkCell ));
746 if ( GetPoints( myCurrElement, P ))
747 aVal = Round( GetValue( P ));
752 double AspectRatio::GetValue( const TSequenceOfXYZ& P )
754 // According to "Mesh quality control" by Nadir Bouhamau referring to
755 // Pascal Jean Frey and Paul-Louis George. Maillages, applications aux elements finis.
756 // Hermes Science publications, Paris 1999 ISBN 2-7462-0024-4
759 int nbNodes = P.size();
764 // Compute aspect ratio
766 if ( nbNodes == 3 ) {
767 // Compute lengths of the sides
768 std::vector< double > aLen (nbNodes);
769 for ( int i = 0; i < nbNodes - 1; i++ )
770 aLen[ i ] = getDistance( P( i + 1 ), P( i + 2 ) );
771 aLen[ nbNodes - 1 ] = getDistance( P( 1 ), P( nbNodes ) );
772 // Q = alfa * h * p / S, where
774 // alfa = sqrt( 3 ) / 6
775 // h - length of the longest edge
776 // p - half perimeter
777 // S - triangle surface
778 const double alfa = sqrt( 3. ) / 6.;
779 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
780 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
781 double anArea = getArea( P( 1 ), P( 2 ), P( 3 ) );
782 if ( anArea <= theEps )
784 return alfa * maxLen * half_perimeter / anArea;
786 else if ( nbNodes == 6 ) { // quadratic triangles
787 // Compute lengths of the sides
788 std::vector< double > aLen (3);
789 aLen[0] = getDistance( P(1), P(3) );
790 aLen[1] = getDistance( P(3), P(5) );
791 aLen[2] = getDistance( P(5), P(1) );
792 // Q = alfa * h * p / S, where
794 // alfa = sqrt( 3 ) / 6
795 // h - length of the longest edge
796 // p - half perimeter
797 // S - triangle surface
798 const double alfa = sqrt( 3. ) / 6.;
799 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
800 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
801 double anArea = getArea( P(1), P(3), P(5) );
802 if ( anArea <= theEps )
804 return alfa * maxLen * half_perimeter / anArea;
806 else if( nbNodes == 4 ) { // quadrangle
807 // Compute lengths of the sides
808 std::vector< double > aLen (4);
809 aLen[0] = getDistance( P(1), P(2) );
810 aLen[1] = getDistance( P(2), P(3) );
811 aLen[2] = getDistance( P(3), P(4) );
812 aLen[3] = getDistance( P(4), P(1) );
813 // Compute lengths of the diagonals
814 std::vector< double > aDia (2);
815 aDia[0] = getDistance( P(1), P(3) );
816 aDia[1] = getDistance( P(2), P(4) );
817 // Compute areas of all triangles which can be built
818 // taking three nodes of the quadrangle
819 std::vector< double > anArea (4);
820 anArea[0] = getArea( P(1), P(2), P(3) );
821 anArea[1] = getArea( P(1), P(2), P(4) );
822 anArea[2] = getArea( P(1), P(3), P(4) );
823 anArea[3] = getArea( P(2), P(3), P(4) );
824 // Q = alpha * L * C1 / C2, where
826 // alpha = sqrt( 1/32 )
827 // L = max( L1, L2, L3, L4, D1, D2 )
828 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
829 // C2 = min( S1, S2, S3, S4 )
830 // Li - lengths of the edges
831 // Di - lengths of the diagonals
832 // Si - areas of the triangles
833 const double alpha = sqrt( 1 / 32. );
834 double L = Max( aLen[ 0 ],
838 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
839 double C1 = sqrt( ( aLen[0] * aLen[0] +
842 aLen[3] * aLen[3] ) / 4. );
843 double C2 = Min( anArea[ 0 ],
845 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
848 return alpha * L * C1 / C2;
850 else if( nbNodes == 8 || nbNodes == 9 ) { // nbNodes==8 - quadratic quadrangle
851 // Compute lengths of the sides
852 std::vector< double > aLen (4);
853 aLen[0] = getDistance( P(1), P(3) );
854 aLen[1] = getDistance( P(3), P(5) );
855 aLen[2] = getDistance( P(5), P(7) );
856 aLen[3] = getDistance( P(7), P(1) );
857 // Compute lengths of the diagonals
858 std::vector< double > aDia (2);
859 aDia[0] = getDistance( P(1), P(5) );
860 aDia[1] = getDistance( P(3), P(7) );
861 // Compute areas of all triangles which can be built
862 // taking three nodes of the quadrangle
863 std::vector< double > anArea (4);
864 anArea[0] = getArea( P(1), P(3), P(5) );
865 anArea[1] = getArea( P(1), P(3), P(7) );
866 anArea[2] = getArea( P(1), P(5), P(7) );
867 anArea[3] = getArea( P(3), P(5), P(7) );
868 // Q = alpha * L * C1 / C2, where
870 // alpha = sqrt( 1/32 )
871 // L = max( L1, L2, L3, L4, D1, D2 )
872 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
873 // C2 = min( S1, S2, S3, S4 )
874 // Li - lengths of the edges
875 // Di - lengths of the diagonals
876 // Si - areas of the triangles
877 const double alpha = sqrt( 1 / 32. );
878 double L = Max( aLen[ 0 ],
882 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
883 double C1 = sqrt( ( aLen[0] * aLen[0] +
886 aLen[3] * aLen[3] ) / 4. );
887 double C2 = Min( anArea[ 0 ],
889 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
892 return alpha * L * C1 / C2;
897 double AspectRatio::GetBadRate( double Value, int /*nbNodes*/ ) const
899 // the aspect ratio is in the range [1.0,infinity]
900 // < 1.0 = very bad, zero area
903 return ( Value < 0.9 ) ? 1000 : Value / 1000.;
906 SMDSAbs_ElementType AspectRatio::GetType() const
912 //================================================================================
914 Class : AspectRatio3D
915 Description : Functor for calculating aspect ratio
917 //================================================================================
921 inline double getHalfPerimeter(double theTria[3]){
922 return (theTria[0] + theTria[1] + theTria[2])/2.0;
925 inline double getArea(double theHalfPerim, double theTria[3]){
926 return sqrt(theHalfPerim*
927 (theHalfPerim-theTria[0])*
928 (theHalfPerim-theTria[1])*
929 (theHalfPerim-theTria[2]));
932 inline double getVolume(double theLen[6]){
933 double a2 = theLen[0]*theLen[0];
934 double b2 = theLen[1]*theLen[1];
935 double c2 = theLen[2]*theLen[2];
936 double d2 = theLen[3]*theLen[3];
937 double e2 = theLen[4]*theLen[4];
938 double f2 = theLen[5]*theLen[5];
939 double P = 4.0*a2*b2*d2;
940 double Q = a2*(b2+d2-e2)-b2*(a2+d2-f2)-d2*(a2+b2-c2);
941 double R = (b2+d2-e2)*(a2+d2-f2)*(a2+d2-f2);
942 return sqrt(P-Q+R)/12.0;
945 inline double getVolume2(double theLen[6]){
946 double a2 = theLen[0]*theLen[0];
947 double b2 = theLen[1]*theLen[1];
948 double c2 = theLen[2]*theLen[2];
949 double d2 = theLen[3]*theLen[3];
950 double e2 = theLen[4]*theLen[4];
951 double f2 = theLen[5]*theLen[5];
953 double P = a2*e2*(b2+c2+d2+f2-a2-e2);
954 double Q = b2*f2*(a2+c2+d2+e2-b2-f2);
955 double R = c2*d2*(a2+b2+e2+f2-c2-d2);
956 double S = a2*b2*d2+b2*c2*e2+a2*c2*f2+d2*e2*f2;
958 return sqrt(P+Q+R-S)/12.0;
961 inline double getVolume(const TSequenceOfXYZ& P){
962 gp_Vec aVec1( P( 2 ) - P( 1 ) );
963 gp_Vec aVec2( P( 3 ) - P( 1 ) );
964 gp_Vec aVec3( P( 4 ) - P( 1 ) );
965 gp_Vec anAreaVec( aVec1 ^ aVec2 );
966 return fabs(aVec3 * anAreaVec) / 6.0;
969 inline double getMaxHeight(double theLen[6])
971 double aHeight = std::max(theLen[0],theLen[1]);
972 aHeight = std::max(aHeight,theLen[2]);
973 aHeight = std::max(aHeight,theLen[3]);
974 aHeight = std::max(aHeight,theLen[4]);
975 aHeight = std::max(aHeight,theLen[5]);
981 double AspectRatio3D::GetValue( long theId )
984 myCurrElement = myMesh->FindElement( theId );
985 if ( myCurrElement && myCurrElement->GetVtkType() == VTK_TETRA )
987 // Action from CoTech | ACTION 31.3:
988 // EURIWARE BO: Homogenize the formulas used to calculate the Controls in SMESH to fit with
989 // those of ParaView. The library used by ParaView for those calculations can be reused in SMESH.
990 vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myCurrElement->getMeshId()]->getGrid();
991 if ( vtkCell* avtkCell = grid->GetCell( myCurrElement->getVtkId() ))
992 aVal = Round( vtkMeshQuality::TetAspectRatio( avtkCell ));
997 if ( GetPoints( myCurrElement, P ))
998 aVal = Round( GetValue( P ));
1003 double AspectRatio3D::GetValue( const TSequenceOfXYZ& P )
1005 double aQuality = 0.0;
1006 if(myCurrElement->IsPoly()) return aQuality;
1008 int nbNodes = P.size();
1010 if(myCurrElement->IsQuadratic()) {
1011 if(nbNodes==10) nbNodes=4; // quadratic tetrahedron
1012 else if(nbNodes==13) nbNodes=5; // quadratic pyramid
1013 else if(nbNodes==15) nbNodes=6; // quadratic pentahedron
1014 else if(nbNodes==20) nbNodes=8; // quadratic hexahedron
1015 else if(nbNodes==27) nbNodes=8; // quadratic hexahedron
1016 else return aQuality;
1022 getDistance(P( 1 ),P( 2 )), // a
1023 getDistance(P( 2 ),P( 3 )), // b
1024 getDistance(P( 3 ),P( 1 )), // c
1025 getDistance(P( 2 ),P( 4 )), // d
1026 getDistance(P( 3 ),P( 4 )), // e
1027 getDistance(P( 1 ),P( 4 )) // f
1029 double aTria[4][3] = {
1030 {aLen[0],aLen[1],aLen[2]}, // abc
1031 {aLen[0],aLen[3],aLen[5]}, // adf
1032 {aLen[1],aLen[3],aLen[4]}, // bde
1033 {aLen[2],aLen[4],aLen[5]} // cef
1035 double aSumArea = 0.0;
1036 double aHalfPerimeter = getHalfPerimeter(aTria[0]);
1037 double anArea = getArea(aHalfPerimeter,aTria[0]);
1039 aHalfPerimeter = getHalfPerimeter(aTria[1]);
1040 anArea = getArea(aHalfPerimeter,aTria[1]);
1042 aHalfPerimeter = getHalfPerimeter(aTria[2]);
1043 anArea = getArea(aHalfPerimeter,aTria[2]);
1045 aHalfPerimeter = getHalfPerimeter(aTria[3]);
1046 anArea = getArea(aHalfPerimeter,aTria[3]);
1048 double aVolume = getVolume(P);
1049 //double aVolume = getVolume(aLen);
1050 double aHeight = getMaxHeight(aLen);
1051 static double aCoeff = sqrt(2.0)/12.0;
1052 if ( aVolume > DBL_MIN )
1053 aQuality = aCoeff*aHeight*aSumArea/aVolume;
1058 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 3 ),P( 5 )};
1059 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1062 gp_XYZ aXYZ[4] = {P( 1 ),P( 3 ),P( 4 ),P( 5 )};
1063 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1066 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 5 )};
1067 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1070 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 4 ),P( 5 )};
1071 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1077 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 6 )};
1078 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1081 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 3 )};
1082 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1085 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 6 )};
1086 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1089 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1090 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1093 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 6 )};
1094 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1097 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 3 )};
1098 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1104 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1105 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1108 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 4 )};
1109 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1112 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 7 )};
1113 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1116 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 8 )};
1117 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1120 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 3 )};
1121 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1124 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 4 )};
1125 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1128 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 7 )};
1129 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1132 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 8 )};
1133 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1136 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 3 )};
1137 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1140 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 4 )};
1141 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1144 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 7 )};
1145 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1148 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 8 )};
1149 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1152 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 1 )};
1153 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1156 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 2 )};
1157 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1160 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 5 )};
1161 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1164 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 6 )};
1165 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1168 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 1 )};
1169 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1172 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 2 )};
1173 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1176 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 5 )};
1177 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1180 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 6 )};
1181 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1184 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 1 )};
1185 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1188 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1189 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1192 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 5 )};
1193 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1196 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 6 )};
1197 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1200 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1201 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1204 gp_XYZ aXYZ[4] = {P( 4 ),P( 5 ),P( 8 ),P( 2 )};
1205 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1208 gp_XYZ aXYZ[4] = {P( 1 ),P( 4 ),P( 5 ),P( 3 )};
1209 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1212 gp_XYZ aXYZ[4] = {P( 3 ),P( 6 ),P( 7 ),P( 1 )};
1213 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1216 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 6 ),P( 4 )};
1217 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1220 gp_XYZ aXYZ[4] = {P( 5 ),P( 6 ),P( 8 ),P( 3 )};
1221 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1224 gp_XYZ aXYZ[4] = {P( 7 ),P( 8 ),P( 6 ),P( 1 )};
1225 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1228 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 7 )};
1229 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1232 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 2 ),P( 5 )};
1233 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1239 gp_XYZ aXYZ[8] = {P( 1 ),P( 2 ),P( 4 ),P( 5 ),P( 7 ),P( 8 ),P( 10 ),P( 11 )};
1240 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1243 gp_XYZ aXYZ[8] = {P( 2 ),P( 3 ),P( 5 ),P( 6 ),P( 8 ),P( 9 ),P( 11 ),P( 12 )};
1244 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1247 gp_XYZ aXYZ[8] = {P( 3 ),P( 4 ),P( 6 ),P( 1 ),P( 9 ),P( 10 ),P( 12 ),P( 7 )};
1248 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1251 } // switch(nbNodes)
1253 if ( nbNodes > 4 ) {
1254 // avaluate aspect ratio of quadranle faces
1255 AspectRatio aspect2D;
1256 SMDS_VolumeTool::VolumeType type = SMDS_VolumeTool::GetType( nbNodes );
1257 int nbFaces = SMDS_VolumeTool::NbFaces( type );
1258 TSequenceOfXYZ points(4);
1259 for ( int i = 0; i < nbFaces; ++i ) { // loop on faces of a volume
1260 if ( SMDS_VolumeTool::NbFaceNodes( type, i ) != 4 )
1262 const int* pInd = SMDS_VolumeTool::GetFaceNodesIndices( type, i, true );
1263 for ( int p = 0; p < 4; ++p ) // loop on nodes of a quadranle face
1264 points( p + 1 ) = P( pInd[ p ] + 1 );
1265 aQuality = std::max( aQuality, aspect2D.GetValue( points ));
1271 double AspectRatio3D::GetBadRate( double Value, int /*nbNodes*/ ) const
1273 // the aspect ratio is in the range [1.0,infinity]
1276 return Value / 1000.;
1279 SMDSAbs_ElementType AspectRatio3D::GetType() const
1281 return SMDSAbs_Volume;
1285 //================================================================================
1288 Description : Functor for calculating warping
1290 //================================================================================
1292 double Warping::GetValue( const TSequenceOfXYZ& P )
1294 if ( P.size() != 4 )
1297 gp_XYZ G = ( P( 1 ) + P( 2 ) + P( 3 ) + P( 4 ) ) / 4.;
1299 double A1 = ComputeA( P( 1 ), P( 2 ), P( 3 ), G );
1300 double A2 = ComputeA( P( 2 ), P( 3 ), P( 4 ), G );
1301 double A3 = ComputeA( P( 3 ), P( 4 ), P( 1 ), G );
1302 double A4 = ComputeA( P( 4 ), P( 1 ), P( 2 ), G );
1304 double val = Max( Max( A1, A2 ), Max( A3, A4 ) );
1306 const double eps = 0.1; // val is in degrees
1308 return val < eps ? 0. : val;
1311 double Warping::ComputeA( const gp_XYZ& thePnt1,
1312 const gp_XYZ& thePnt2,
1313 const gp_XYZ& thePnt3,
1314 const gp_XYZ& theG ) const
1316 double aLen1 = gp_Pnt( thePnt1 ).Distance( gp_Pnt( thePnt2 ) );
1317 double aLen2 = gp_Pnt( thePnt2 ).Distance( gp_Pnt( thePnt3 ) );
1318 double L = Min( aLen1, aLen2 ) * 0.5;
1322 gp_XYZ GI = ( thePnt2 + thePnt1 ) / 2. - theG;
1323 gp_XYZ GJ = ( thePnt3 + thePnt2 ) / 2. - theG;
1324 gp_XYZ N = GI.Crossed( GJ );
1326 if ( N.Modulus() < gp::Resolution() )
1331 double H = ( thePnt2 - theG ).Dot( N );
1332 return asin( fabs( H / L ) ) * 180. / M_PI;
1335 double Warping::GetBadRate( double Value, int /*nbNodes*/ ) const
1337 // the warp is in the range [0.0,PI/2]
1338 // 0.0 = good (no warp)
1339 // PI/2 = bad (face pliee)
1343 SMDSAbs_ElementType Warping::GetType() const
1345 return SMDSAbs_Face;
1349 //================================================================================
1352 Description : Functor for calculating taper
1354 //================================================================================
1356 double Taper::GetValue( const TSequenceOfXYZ& P )
1358 if ( P.size() != 4 )
1362 double J1 = getArea( P( 4 ), P( 1 ), P( 2 ) );
1363 double J2 = getArea( P( 3 ), P( 1 ), P( 2 ) );
1364 double J3 = getArea( P( 2 ), P( 3 ), P( 4 ) );
1365 double J4 = getArea( P( 3 ), P( 4 ), P( 1 ) );
1367 double JA = 0.25 * ( J1 + J2 + J3 + J4 );
1371 double T1 = fabs( ( J1 - JA ) / JA );
1372 double T2 = fabs( ( J2 - JA ) / JA );
1373 double T3 = fabs( ( J3 - JA ) / JA );
1374 double T4 = fabs( ( J4 - JA ) / JA );
1376 double val = Max( Max( T1, T2 ), Max( T3, T4 ) );
1378 const double eps = 0.01;
1380 return val < eps ? 0. : val;
1383 double Taper::GetBadRate( double Value, int /*nbNodes*/ ) const
1385 // the taper is in the range [0.0,1.0]
1386 // 0.0 = good (no taper)
1387 // 1.0 = bad (les cotes opposes sont allignes)
1391 SMDSAbs_ElementType Taper::GetType() const
1393 return SMDSAbs_Face;
1396 //================================================================================
1399 Description : Functor for calculating skew in degrees
1401 //================================================================================
1403 static inline double skewAngle( const gp_XYZ& p1, const gp_XYZ& p2, const gp_XYZ& p3 )
1405 gp_XYZ p12 = ( p2 + p1 ) / 2.;
1406 gp_XYZ p23 = ( p3 + p2 ) / 2.;
1407 gp_XYZ p31 = ( p3 + p1 ) / 2.;
1409 gp_Vec v1( p31 - p2 ), v2( p12 - p23 );
1411 return v1.Magnitude() < gp::Resolution() || v2.Magnitude() < gp::Resolution() ? 0. : v1.Angle( v2 );
1414 double Skew::GetValue( const TSequenceOfXYZ& P )
1416 if ( P.size() != 3 && P.size() != 4 )
1420 const double PI2 = M_PI / 2.;
1421 if ( P.size() == 3 )
1423 double A0 = fabs( PI2 - skewAngle( P( 3 ), P( 1 ), P( 2 ) ) );
1424 double A1 = fabs( PI2 - skewAngle( P( 1 ), P( 2 ), P( 3 ) ) );
1425 double A2 = fabs( PI2 - skewAngle( P( 2 ), P( 3 ), P( 1 ) ) );
1427 return Max( A0, Max( A1, A2 ) ) * 180. / M_PI;
1431 gp_XYZ p12 = ( P( 1 ) + P( 2 ) ) / 2.;
1432 gp_XYZ p23 = ( P( 2 ) + P( 3 ) ) / 2.;
1433 gp_XYZ p34 = ( P( 3 ) + P( 4 ) ) / 2.;
1434 gp_XYZ p41 = ( P( 4 ) + P( 1 ) ) / 2.;
1436 gp_Vec v1( p34 - p12 ), v2( p23 - p41 );
1437 double A = v1.Magnitude() <= gp::Resolution() || v2.Magnitude() <= gp::Resolution()
1438 ? 0. : fabs( PI2 - v1.Angle( v2 ) );
1440 double val = A * 180. / M_PI;
1442 const double eps = 0.1; // val is in degrees
1444 return val < eps ? 0. : val;
1448 double Skew::GetBadRate( double Value, int /*nbNodes*/ ) const
1450 // the skew is in the range [0.0,PI/2].
1456 SMDSAbs_ElementType Skew::GetType() const
1458 return SMDSAbs_Face;
1462 //================================================================================
1465 Description : Functor for calculating area
1467 //================================================================================
1469 double Area::GetValue( const TSequenceOfXYZ& P )
1474 gp_Vec aVec1( P(2) - P(1) );
1475 gp_Vec aVec2( P(3) - P(1) );
1476 gp_Vec SumVec = aVec1 ^ aVec2;
1478 for (int i=4; i<=P.size(); i++)
1480 gp_Vec aVec1( P(i-1) - P(1) );
1481 gp_Vec aVec2( P(i) - P(1) );
1482 gp_Vec tmp = aVec1 ^ aVec2;
1485 val = SumVec.Magnitude() * 0.5;
1490 double Area::GetBadRate( double Value, int /*nbNodes*/ ) const
1492 // meaningless as it is not a quality control functor
1496 SMDSAbs_ElementType Area::GetType() const
1498 return SMDSAbs_Face;
1501 //================================================================================
1504 Description : Functor for calculating length of edge
1506 //================================================================================
1508 double Length::GetValue( const TSequenceOfXYZ& P )
1510 switch ( P.size() ) {
1511 case 2: return getDistance( P( 1 ), P( 2 ) );
1512 case 3: return getDistance( P( 1 ), P( 2 ) ) + getDistance( P( 2 ), P( 3 ) );
1517 double Length::GetBadRate( double Value, int /*nbNodes*/ ) const
1519 // meaningless as it is not quality control functor
1523 SMDSAbs_ElementType Length::GetType() const
1525 return SMDSAbs_Edge;
1528 //================================================================================
1531 Description : Functor for calculating minimal length of edge
1533 //================================================================================
1535 double Length2D::GetValue( long theElementId )
1539 if ( GetPoints( theElementId, P ))
1543 SMDSAbs_EntityType aType = P.getElementEntity();
1546 case SMDSEntity_Edge:
1548 aVal = getDistance( P( 1 ), P( 2 ) );
1550 case SMDSEntity_Quad_Edge:
1551 if (len == 3) // quadratic edge
1552 aVal = getDistance(P( 1 ),P( 3 )) + getDistance(P( 3 ),P( 2 ));
1554 case SMDSEntity_Triangle:
1555 if (len == 3){ // triangles
1556 double L1 = getDistance(P( 1 ),P( 2 ));
1557 double L2 = getDistance(P( 2 ),P( 3 ));
1558 double L3 = getDistance(P( 3 ),P( 1 ));
1559 aVal = Min(L1,Min(L2,L3));
1562 case SMDSEntity_Quadrangle:
1563 if (len == 4){ // quadrangles
1564 double L1 = getDistance(P( 1 ),P( 2 ));
1565 double L2 = getDistance(P( 2 ),P( 3 ));
1566 double L3 = getDistance(P( 3 ),P( 4 ));
1567 double L4 = getDistance(P( 4 ),P( 1 ));
1568 aVal = Min(Min(L1,L2),Min(L3,L4));
1571 case SMDSEntity_Quad_Triangle:
1572 case SMDSEntity_BiQuad_Triangle:
1573 if (len >= 6){ // quadratic triangles
1574 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1575 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1576 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
1577 aVal = Min(L1,Min(L2,L3));
1580 case SMDSEntity_Quad_Quadrangle:
1581 case SMDSEntity_BiQuad_Quadrangle:
1582 if (len >= 8){ // quadratic quadrangles
1583 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1584 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1585 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
1586 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
1587 aVal = Min(Min(L1,L2),Min(L3,L4));
1590 case SMDSEntity_Tetra:
1591 if (len == 4){ // tetrahedra
1592 double L1 = getDistance(P( 1 ),P( 2 ));
1593 double L2 = getDistance(P( 2 ),P( 3 ));
1594 double L3 = getDistance(P( 3 ),P( 1 ));
1595 double L4 = getDistance(P( 1 ),P( 4 ));
1596 double L5 = getDistance(P( 2 ),P( 4 ));
1597 double L6 = getDistance(P( 3 ),P( 4 ));
1598 aVal = Min(Min(Min(L1,L2),Min(L3,L4)),Min(L5,L6));
1601 case SMDSEntity_Pyramid:
1602 if (len == 5){ // piramids
1603 double L1 = getDistance(P( 1 ),P( 2 ));
1604 double L2 = getDistance(P( 2 ),P( 3 ));
1605 double L3 = getDistance(P( 3 ),P( 4 ));
1606 double L4 = getDistance(P( 4 ),P( 1 ));
1607 double L5 = getDistance(P( 1 ),P( 5 ));
1608 double L6 = getDistance(P( 2 ),P( 5 ));
1609 double L7 = getDistance(P( 3 ),P( 5 ));
1610 double L8 = getDistance(P( 4 ),P( 5 ));
1612 aVal = Min(Min(Min(L1,L2),Min(L3,L4)),Min(L5,L6));
1613 aVal = Min(aVal,Min(L7,L8));
1616 case SMDSEntity_Penta:
1617 if (len == 6) { // pentaidres
1618 double L1 = getDistance(P( 1 ),P( 2 ));
1619 double L2 = getDistance(P( 2 ),P( 3 ));
1620 double L3 = getDistance(P( 3 ),P( 1 ));
1621 double L4 = getDistance(P( 4 ),P( 5 ));
1622 double L5 = getDistance(P( 5 ),P( 6 ));
1623 double L6 = getDistance(P( 6 ),P( 4 ));
1624 double L7 = getDistance(P( 1 ),P( 4 ));
1625 double L8 = getDistance(P( 2 ),P( 5 ));
1626 double L9 = getDistance(P( 3 ),P( 6 ));
1628 aVal = Min(Min(Min(L1,L2),Min(L3,L4)),Min(L5,L6));
1629 aVal = Min(aVal,Min(Min(L7,L8),L9));
1632 case SMDSEntity_Hexa:
1633 if (len == 8){ // hexahedron
1634 double L1 = getDistance(P( 1 ),P( 2 ));
1635 double L2 = getDistance(P( 2 ),P( 3 ));
1636 double L3 = getDistance(P( 3 ),P( 4 ));
1637 double L4 = getDistance(P( 4 ),P( 1 ));
1638 double L5 = getDistance(P( 5 ),P( 6 ));
1639 double L6 = getDistance(P( 6 ),P( 7 ));
1640 double L7 = getDistance(P( 7 ),P( 8 ));
1641 double L8 = getDistance(P( 8 ),P( 5 ));
1642 double L9 = getDistance(P( 1 ),P( 5 ));
1643 double L10= getDistance(P( 2 ),P( 6 ));
1644 double L11= getDistance(P( 3 ),P( 7 ));
1645 double L12= getDistance(P( 4 ),P( 8 ));
1647 aVal = Min(Min(Min(L1,L2),Min(L3,L4)),Min(L5,L6));
1648 aVal = Min(aVal,Min(Min(L7,L8),Min(L9,L10)));
1649 aVal = Min(aVal,Min(L11,L12));
1652 case SMDSEntity_Quad_Tetra:
1653 if (len == 10){ // quadratic tetraidrs
1654 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
1655 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
1656 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
1657 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1658 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
1659 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
1660 aVal = Min(Min(Min(L1,L2),Min(L3,L4)),Min(L5,L6));
1663 case SMDSEntity_Quad_Pyramid:
1664 if (len == 13){ // quadratic piramids
1665 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
1666 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
1667 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1668 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1669 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1670 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
1671 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
1672 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
1673 aVal = Min(Min(Min(L1,L2),Min(L3,L4)),Min(L5,L6));
1674 aVal = Min(aVal,Min(L7,L8));
1677 case SMDSEntity_Quad_Penta:
1678 if (len == 15){ // quadratic pentaidres
1679 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
1680 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
1681 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1682 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1683 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
1684 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
1685 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
1686 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
1687 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
1688 aVal = Min(Min(Min(L1,L2),Min(L3,L4)),Min(L5,L6));
1689 aVal = Min(aVal,Min(Min(L7,L8),L9));
1692 case SMDSEntity_Quad_Hexa:
1693 case SMDSEntity_TriQuad_Hexa:
1694 if (len >= 20) { // quadratic hexaider
1695 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
1696 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
1697 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
1698 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
1699 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
1700 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
1701 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
1702 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
1703 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
1704 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
1705 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
1706 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
1707 aVal = Min(Min(Min(L1,L2),Min(L3,L4)),Min(L5,L6));
1708 aVal = Min(aVal,Min(Min(L7,L8),Min(L9,L10)));
1709 aVal = Min(aVal,Min(L11,L12));
1712 case SMDSEntity_Polygon:
1714 aVal = getDistance( P(1), P( P.size() ));
1715 for ( size_t i = 1; i < P.size(); ++i )
1716 aVal = Min( aVal, getDistance( P( i ), P( i+1 )));
1719 case SMDSEntity_Quad_Polygon:
1721 aVal = getDistance( P(1), P( P.size() )) + getDistance( P(P.size()), P( P.size()-1 ));
1722 for ( size_t i = 1; i < P.size()-1; i += 2 )
1723 aVal = Min( aVal, getDistance( P( i ), P( i+1 )) + getDistance( P( i+1 ), P( i+2 )));
1726 case SMDSEntity_Hexagonal_Prism:
1727 if (len == 12) { // hexagonal prism
1728 double L1 = getDistance(P( 1 ),P( 2 ));
1729 double L2 = getDistance(P( 2 ),P( 3 ));
1730 double L3 = getDistance(P( 3 ),P( 4 ));
1731 double L4 = getDistance(P( 4 ),P( 5 ));
1732 double L5 = getDistance(P( 5 ),P( 6 ));
1733 double L6 = getDistance(P( 6 ),P( 1 ));
1735 double L7 = getDistance(P( 7 ), P( 8 ));
1736 double L8 = getDistance(P( 8 ), P( 9 ));
1737 double L9 = getDistance(P( 9 ), P( 10 ));
1738 double L10= getDistance(P( 10 ),P( 11 ));
1739 double L11= getDistance(P( 11 ),P( 12 ));
1740 double L12= getDistance(P( 12 ),P( 7 ));
1742 double L13 = getDistance(P( 1 ),P( 7 ));
1743 double L14 = getDistance(P( 2 ),P( 8 ));
1744 double L15 = getDistance(P( 3 ),P( 9 ));
1745 double L16 = getDistance(P( 4 ),P( 10 ));
1746 double L17 = getDistance(P( 5 ),P( 11 ));
1747 double L18 = getDistance(P( 6 ),P( 12 ));
1748 aVal = Min(Min(Min(L1,L2),Min(L3,L4)),Min(L5,L6));
1749 aVal = Min(aVal, Min(Min(Min(L7,L8),Min(L9,L10)),Min(L11,L12)));
1750 aVal = Min(aVal, Min(Min(Min(L13,L14),Min(L15,L16)),Min(L17,L18)));
1753 case SMDSEntity_Polyhedra:
1765 if ( myPrecision >= 0 )
1767 double prec = pow( 10., (double)( myPrecision ) );
1768 aVal = floor( aVal * prec + 0.5 ) / prec;
1777 double Length2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1779 // meaningless as it is not a quality control functor
1783 SMDSAbs_ElementType Length2D::GetType() const
1785 return SMDSAbs_Face;
1788 Length2D::Value::Value(double theLength,long thePntId1, long thePntId2):
1791 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1792 if(thePntId1 > thePntId2){
1793 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1797 bool Length2D::Value::operator<(const Length2D::Value& x) const
1799 if(myPntId[0] < x.myPntId[0]) return true;
1800 if(myPntId[0] == x.myPntId[0])
1801 if(myPntId[1] < x.myPntId[1]) return true;
1805 void Length2D::GetValues(TValues& theValues)
1808 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1809 for(; anIter->more(); ){
1810 const SMDS_MeshFace* anElem = anIter->next();
1812 if(anElem->IsQuadratic()) {
1813 const SMDS_VtkFace* F =
1814 dynamic_cast<const SMDS_VtkFace*>(anElem);
1815 // use special nodes iterator
1816 SMDS_ElemIteratorPtr anIter = F->interlacedNodesElemIterator();
1821 const SMDS_MeshElement* aNode;
1823 aNode = anIter->next();
1824 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1825 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1826 aNodeId[0] = aNodeId[1] = aNode->GetID();
1829 for(; anIter->more(); ){
1830 const SMDS_MeshNode* N1 = static_cast<const SMDS_MeshNode*> (anIter->next());
1831 P[2] = gp_Pnt(N1->X(),N1->Y(),N1->Z());
1832 aNodeId[2] = N1->GetID();
1833 aLength = P[1].Distance(P[2]);
1834 if(!anIter->more()) break;
1835 const SMDS_MeshNode* N2 = static_cast<const SMDS_MeshNode*> (anIter->next());
1836 P[3] = gp_Pnt(N2->X(),N2->Y(),N2->Z());
1837 aNodeId[3] = N2->GetID();
1838 aLength += P[2].Distance(P[3]);
1839 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1840 Value aValue2(aLength,aNodeId[2],aNodeId[3]);
1842 aNodeId[1] = aNodeId[3];
1843 theValues.insert(aValue1);
1844 theValues.insert(aValue2);
1846 aLength += P[2].Distance(P[0]);
1847 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1848 Value aValue2(aLength,aNodeId[2],aNodeId[0]);
1849 theValues.insert(aValue1);
1850 theValues.insert(aValue2);
1853 SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
1858 const SMDS_MeshElement* aNode;
1859 if(aNodesIter->more()){
1860 aNode = aNodesIter->next();
1861 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1862 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1863 aNodeId[0] = aNodeId[1] = aNode->GetID();
1866 for(; aNodesIter->more(); ){
1867 aNode = aNodesIter->next();
1868 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1869 long anId = aNode->GetID();
1871 P[2] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1873 aLength = P[1].Distance(P[2]);
1875 Value aValue(aLength,aNodeId[1],anId);
1878 theValues.insert(aValue);
1881 aLength = P[0].Distance(P[1]);
1883 Value aValue(aLength,aNodeId[0],aNodeId[1]);
1884 theValues.insert(aValue);
1889 //================================================================================
1891 Class : MultiConnection
1892 Description : Functor for calculating number of faces conneted to the edge
1894 //================================================================================
1896 double MultiConnection::GetValue( const TSequenceOfXYZ& P )
1900 double MultiConnection::GetValue( long theId )
1902 return getNbMultiConnection( myMesh, theId );
1905 double MultiConnection::GetBadRate( double Value, int /*nbNodes*/ ) const
1907 // meaningless as it is not quality control functor
1911 SMDSAbs_ElementType MultiConnection::GetType() const
1913 return SMDSAbs_Edge;
1916 //================================================================================
1918 Class : MultiConnection2D
1919 Description : Functor for calculating number of faces conneted to the edge
1921 //================================================================================
1923 double MultiConnection2D::GetValue( const TSequenceOfXYZ& P )
1928 double MultiConnection2D::GetValue( long theElementId )
1932 const SMDS_MeshElement* aFaceElem = myMesh->FindElement(theElementId);
1933 SMDSAbs_ElementType aType = aFaceElem->GetType();
1938 int i = 0, len = aFaceElem->NbNodes();
1939 SMDS_ElemIteratorPtr anIter = aFaceElem->nodesIterator();
1942 const SMDS_MeshNode *aNode, *aNode0;
1943 TColStd_MapOfInteger aMap, aMapPrev;
1945 for (i = 0; i <= len; i++) {
1950 if (anIter->more()) {
1951 aNode = (SMDS_MeshNode*)anIter->next();
1959 if (i == 0) aNode0 = aNode;
1961 SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
1962 while (anElemIter->more()) {
1963 const SMDS_MeshElement* anElem = anElemIter->next();
1964 if (anElem != 0 && anElem->GetType() == SMDSAbs_Face) {
1965 int anId = anElem->GetID();
1968 if (aMapPrev.Contains(anId)) {
1973 aResult = Max(aResult, aNb);
1984 double MultiConnection2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1986 // meaningless as it is not quality control functor
1990 SMDSAbs_ElementType MultiConnection2D::GetType() const
1992 return SMDSAbs_Face;
1995 MultiConnection2D::Value::Value(long thePntId1, long thePntId2)
1997 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1998 if(thePntId1 > thePntId2){
1999 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
2003 bool MultiConnection2D::Value::operator<(const MultiConnection2D::Value& x) const
2005 if(myPntId[0] < x.myPntId[0]) return true;
2006 if(myPntId[0] == x.myPntId[0])
2007 if(myPntId[1] < x.myPntId[1]) return true;
2011 void MultiConnection2D::GetValues(MValues& theValues)
2013 if ( !myMesh ) return;
2014 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
2015 for(; anIter->more(); ){
2016 const SMDS_MeshFace* anElem = anIter->next();
2017 SMDS_ElemIteratorPtr aNodesIter;
2018 if ( anElem->IsQuadratic() )
2019 aNodesIter = dynamic_cast<const SMDS_VtkFace*>
2020 (anElem)->interlacedNodesElemIterator();
2022 aNodesIter = anElem->nodesIterator();
2025 //int aNbConnects=0;
2026 const SMDS_MeshNode* aNode0;
2027 const SMDS_MeshNode* aNode1;
2028 const SMDS_MeshNode* aNode2;
2029 if(aNodesIter->more()){
2030 aNode0 = (SMDS_MeshNode*) aNodesIter->next();
2032 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode1;
2033 aNodeId[0] = aNodeId[1] = aNodes->GetID();
2035 for(; aNodesIter->more(); ) {
2036 aNode2 = (SMDS_MeshNode*) aNodesIter->next();
2037 long anId = aNode2->GetID();
2040 Value aValue(aNodeId[1],aNodeId[2]);
2041 MValues::iterator aItr = theValues.find(aValue);
2042 if (aItr != theValues.end()){
2047 theValues[aValue] = 1;
2050 //cout << "NodeIds: "<<aNodeId[1]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
2051 aNodeId[1] = aNodeId[2];
2054 Value aValue(aNodeId[0],aNodeId[2]);
2055 MValues::iterator aItr = theValues.find(aValue);
2056 if (aItr != theValues.end()) {
2061 theValues[aValue] = 1;
2064 //cout << "NodeIds: "<<aNodeId[0]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
2069 //================================================================================
2071 Class : BallDiameter
2072 Description : Functor returning diameter of a ball element
2074 //================================================================================
2076 double BallDiameter::GetValue( long theId )
2078 double diameter = 0;
2080 if ( const SMDS_BallElement* ball =
2081 dynamic_cast<const SMDS_BallElement*>( myMesh->FindElement( theId )))
2083 diameter = ball->GetDiameter();
2088 double BallDiameter::GetBadRate( double Value, int /*nbNodes*/ ) const
2090 // meaningless as it is not a quality control functor
2094 SMDSAbs_ElementType BallDiameter::GetType() const
2096 return SMDSAbs_Ball;
2104 //================================================================================
2106 Class : BadOrientedVolume
2107 Description : Predicate bad oriented volumes
2109 //================================================================================
2111 BadOrientedVolume::BadOrientedVolume()
2116 void BadOrientedVolume::SetMesh( const SMDS_Mesh* theMesh )
2121 bool BadOrientedVolume::IsSatisfy( long theId )
2126 SMDS_VolumeTool vTool( myMesh->FindElement( theId ));
2127 return !vTool.IsForward();
2130 SMDSAbs_ElementType BadOrientedVolume::GetType() const
2132 return SMDSAbs_Volume;
2136 Class : BareBorderVolume
2139 bool BareBorderVolume::IsSatisfy(long theElementId )
2141 SMDS_VolumeTool myTool;
2142 if ( myTool.Set( myMesh->FindElement(theElementId)))
2144 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2145 if ( myTool.IsFreeFace( iF ))
2147 const SMDS_MeshNode** n = myTool.GetFaceNodes(iF);
2148 vector< const SMDS_MeshNode*> nodes( n, n+myTool.NbFaceNodes(iF));
2149 if ( !myMesh->FindElement( nodes, SMDSAbs_Face, /*Nomedium=*/false))
2156 //================================================================================
2158 Class : BareBorderFace
2160 //================================================================================
2162 bool BareBorderFace::IsSatisfy(long theElementId )
2165 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2167 if ( face->GetType() == SMDSAbs_Face )
2169 int nbN = face->NbCornerNodes();
2170 for ( int i = 0; i < nbN && !ok; ++i )
2172 // check if a link is shared by another face
2173 const SMDS_MeshNode* n1 = face->GetNode( i );
2174 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2175 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2176 bool isShared = false;
2177 while ( !isShared && fIt->more() )
2179 const SMDS_MeshElement* f = fIt->next();
2180 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2184 const int iQuad = face->IsQuadratic();
2185 myLinkNodes.resize( 2 + iQuad);
2186 myLinkNodes[0] = n1;
2187 myLinkNodes[1] = n2;
2189 myLinkNodes[2] = face->GetNode( i+nbN );
2190 ok = !myMesh->FindElement( myLinkNodes, SMDSAbs_Edge, /*noMedium=*/false);
2198 //================================================================================
2200 Class : OverConstrainedVolume
2202 //================================================================================
2204 bool OverConstrainedVolume::IsSatisfy(long theElementId )
2206 // An element is over-constrained if it has N-1 free borders where
2207 // N is the number of edges/faces for a 2D/3D element.
2208 SMDS_VolumeTool myTool;
2209 if ( myTool.Set( myMesh->FindElement(theElementId)))
2211 int nbSharedFaces = 0;
2212 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2213 if ( !myTool.IsFreeFace( iF ) && ++nbSharedFaces > 1 )
2215 return ( nbSharedFaces == 1 );
2220 //================================================================================
2222 Class : OverConstrainedFace
2224 //================================================================================
2226 bool OverConstrainedFace::IsSatisfy(long theElementId )
2228 // An element is over-constrained if it has N-1 free borders where
2229 // N is the number of edges/faces for a 2D/3D element.
2230 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2231 if ( face->GetType() == SMDSAbs_Face )
2233 int nbSharedBorders = 0;
2234 int nbN = face->NbCornerNodes();
2235 for ( int i = 0; i < nbN; ++i )
2237 // check if a link is shared by another face
2238 const SMDS_MeshNode* n1 = face->GetNode( i );
2239 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2240 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2241 bool isShared = false;
2242 while ( !isShared && fIt->more() )
2244 const SMDS_MeshElement* f = fIt->next();
2245 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2247 if ( isShared && ++nbSharedBorders > 1 )
2250 return ( nbSharedBorders == 1 );
2255 //================================================================================
2257 Class : CoincidentNodes
2258 Description : Predicate of Coincident nodes
2260 //================================================================================
2262 CoincidentNodes::CoincidentNodes()
2267 bool CoincidentNodes::IsSatisfy( long theElementId )
2269 return myCoincidentIDs.Contains( theElementId );
2272 SMDSAbs_ElementType CoincidentNodes::GetType() const
2274 return SMDSAbs_Node;
2277 void CoincidentNodes::SetMesh( const SMDS_Mesh* theMesh )
2279 myMeshModifTracer.SetMesh( theMesh );
2280 if ( myMeshModifTracer.IsMeshModified() )
2282 TIDSortedNodeSet nodesToCheck;
2283 SMDS_NodeIteratorPtr nIt = theMesh->nodesIterator(/*idInceasingOrder=*/true);
2284 while ( nIt->more() )
2285 nodesToCheck.insert( nodesToCheck.end(), nIt->next() );
2287 list< list< const SMDS_MeshNode*> > nodeGroups;
2288 SMESH_OctreeNode::FindCoincidentNodes ( nodesToCheck, &nodeGroups, myToler );
2290 myCoincidentIDs.Clear();
2291 list< list< const SMDS_MeshNode*> >::iterator groupIt = nodeGroups.begin();
2292 for ( ; groupIt != nodeGroups.end(); ++groupIt )
2294 list< const SMDS_MeshNode*>& coincNodes = *groupIt;
2295 list< const SMDS_MeshNode*>::iterator n = coincNodes.begin();
2296 for ( ; n != coincNodes.end(); ++n )
2297 myCoincidentIDs.Add( (*n)->GetID() );
2302 //================================================================================
2304 Class : CoincidentElements
2305 Description : Predicate of Coincident Elements
2306 Note : This class is suitable only for visualization of Coincident Elements
2308 //================================================================================
2310 CoincidentElements::CoincidentElements()
2315 void CoincidentElements::SetMesh( const SMDS_Mesh* theMesh )
2320 bool CoincidentElements::IsSatisfy( long theElementId )
2322 if ( !myMesh ) return false;
2324 if ( const SMDS_MeshElement* e = myMesh->FindElement( theElementId ))
2326 if ( e->GetType() != GetType() ) return false;
2327 set< const SMDS_MeshNode* > elemNodes( e->begin_nodes(), e->end_nodes() );
2328 const int nbNodes = e->NbNodes();
2329 SMDS_ElemIteratorPtr invIt = (*elemNodes.begin())->GetInverseElementIterator( GetType() );
2330 while ( invIt->more() )
2332 const SMDS_MeshElement* e2 = invIt->next();
2333 if ( e2 == e || e2->NbNodes() != nbNodes ) continue;
2335 bool sameNodes = true;
2336 for ( size_t i = 0; i < elemNodes.size() && sameNodes; ++i )
2337 sameNodes = ( elemNodes.count( e2->GetNode( i )));
2345 SMDSAbs_ElementType CoincidentElements1D::GetType() const
2347 return SMDSAbs_Edge;
2349 SMDSAbs_ElementType CoincidentElements2D::GetType() const
2351 return SMDSAbs_Face;
2353 SMDSAbs_ElementType CoincidentElements3D::GetType() const
2355 return SMDSAbs_Volume;
2359 //================================================================================
2362 Description : Predicate for free borders
2364 //================================================================================
2366 FreeBorders::FreeBorders()
2371 void FreeBorders::SetMesh( const SMDS_Mesh* theMesh )
2376 bool FreeBorders::IsSatisfy( long theId )
2378 return getNbMultiConnection( myMesh, theId ) == 1;
2381 SMDSAbs_ElementType FreeBorders::GetType() const
2383 return SMDSAbs_Edge;
2387 //================================================================================
2390 Description : Predicate for free Edges
2392 //================================================================================
2394 FreeEdges::FreeEdges()
2399 void FreeEdges::SetMesh( const SMDS_Mesh* theMesh )
2404 bool FreeEdges::IsFreeEdge( const SMDS_MeshNode** theNodes, const int theFaceId )
2406 TColStd_MapOfInteger aMap;
2407 for ( int i = 0; i < 2; i++ )
2409 SMDS_ElemIteratorPtr anElemIter = theNodes[ i ]->GetInverseElementIterator(SMDSAbs_Face);
2410 while( anElemIter->more() )
2412 if ( const SMDS_MeshElement* anElem = anElemIter->next())
2414 const int anId = anElem->GetID();
2415 if ( anId != theFaceId && !aMap.Add( anId ))
2423 bool FreeEdges::IsSatisfy( long theId )
2428 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2429 if ( aFace == 0 || aFace->GetType() != SMDSAbs_Face || aFace->NbNodes() < 3 )
2432 SMDS_ElemIteratorPtr anIter;
2433 if ( aFace->IsQuadratic() ) {
2434 anIter = dynamic_cast<const SMDS_VtkFace*>
2435 (aFace)->interlacedNodesElemIterator();
2438 anIter = aFace->nodesIterator();
2443 int i = 0, nbNodes = aFace->NbNodes();
2444 std::vector <const SMDS_MeshNode*> aNodes( nbNodes+1 );
2445 while( anIter->more() )
2447 const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
2450 aNodes[ i++ ] = aNode;
2452 aNodes[ nbNodes ] = aNodes[ 0 ];
2454 for ( i = 0; i < nbNodes; i++ )
2455 if ( IsFreeEdge( &aNodes[ i ], theId ) )
2461 SMDSAbs_ElementType FreeEdges::GetType() const
2463 return SMDSAbs_Face;
2466 FreeEdges::Border::Border(long theElemId, long thePntId1, long thePntId2):
2469 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
2470 if(thePntId1 > thePntId2){
2471 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
2475 bool FreeEdges::Border::operator<(const FreeEdges::Border& x) const{
2476 if(myPntId[0] < x.myPntId[0]) return true;
2477 if(myPntId[0] == x.myPntId[0])
2478 if(myPntId[1] < x.myPntId[1]) return true;
2482 inline void UpdateBorders(const FreeEdges::Border& theBorder,
2483 FreeEdges::TBorders& theRegistry,
2484 FreeEdges::TBorders& theContainer)
2486 if(theRegistry.find(theBorder) == theRegistry.end()){
2487 theRegistry.insert(theBorder);
2488 theContainer.insert(theBorder);
2490 theContainer.erase(theBorder);
2494 void FreeEdges::GetBoreders(TBorders& theBorders)
2497 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
2498 for(; anIter->more(); ){
2499 const SMDS_MeshFace* anElem = anIter->next();
2500 long anElemId = anElem->GetID();
2501 SMDS_ElemIteratorPtr aNodesIter;
2502 if ( anElem->IsQuadratic() )
2503 aNodesIter = static_cast<const SMDS_VtkFace*>(anElem)->
2504 interlacedNodesElemIterator();
2506 aNodesIter = anElem->nodesIterator();
2508 const SMDS_MeshElement* aNode;
2509 if(aNodesIter->more()){
2510 aNode = aNodesIter->next();
2511 aNodeId[0] = aNodeId[1] = aNode->GetID();
2513 for(; aNodesIter->more(); ){
2514 aNode = aNodesIter->next();
2515 long anId = aNode->GetID();
2516 Border aBorder(anElemId,aNodeId[1],anId);
2518 UpdateBorders(aBorder,aRegistry,theBorders);
2520 Border aBorder(anElemId,aNodeId[0],aNodeId[1]);
2521 UpdateBorders(aBorder,aRegistry,theBorders);
2525 //================================================================================
2528 Description : Predicate for free nodes
2530 //================================================================================
2532 FreeNodes::FreeNodes()
2537 void FreeNodes::SetMesh( const SMDS_Mesh* theMesh )
2542 bool FreeNodes::IsSatisfy( long theNodeId )
2544 const SMDS_MeshNode* aNode = myMesh->FindNode( theNodeId );
2548 return (aNode->NbInverseElements() < 1);
2551 SMDSAbs_ElementType FreeNodes::GetType() const
2553 return SMDSAbs_Node;
2557 //================================================================================
2560 Description : Predicate for free faces
2562 //================================================================================
2564 FreeFaces::FreeFaces()
2569 void FreeFaces::SetMesh( const SMDS_Mesh* theMesh )
2574 bool FreeFaces::IsSatisfy( long theId )
2576 if (!myMesh) return false;
2577 // check that faces nodes refers to less than two common volumes
2578 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2579 if ( !aFace || aFace->GetType() != SMDSAbs_Face )
2582 int nbNode = aFace->NbNodes();
2584 // collect volumes to check that number of volumes with count equal nbNode not less than 2
2585 typedef map< SMDS_MeshElement*, int > TMapOfVolume; // map of volume counters
2586 typedef map< SMDS_MeshElement*, int >::iterator TItrMapOfVolume; // iterator
2587 TMapOfVolume mapOfVol;
2589 SMDS_ElemIteratorPtr nodeItr = aFace->nodesIterator();
2590 while ( nodeItr->more() ) {
2591 const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>(nodeItr->next());
2592 if ( !aNode ) continue;
2593 SMDS_ElemIteratorPtr volItr = aNode->GetInverseElementIterator(SMDSAbs_Volume);
2594 while ( volItr->more() ) {
2595 SMDS_MeshElement* aVol = (SMDS_MeshElement*)volItr->next();
2596 TItrMapOfVolume itr = mapOfVol.insert(make_pair(aVol, 0)).first;
2601 TItrMapOfVolume volItr = mapOfVol.begin();
2602 TItrMapOfVolume volEnd = mapOfVol.end();
2603 for ( ; volItr != volEnd; ++volItr )
2604 if ( (*volItr).second >= nbNode )
2606 // face is not free if number of volumes constructed on thier nodes more than one
2610 SMDSAbs_ElementType FreeFaces::GetType() const
2612 return SMDSAbs_Face;
2615 //================================================================================
2617 Class : LinearOrQuadratic
2618 Description : Predicate to verify whether a mesh element is linear
2620 //================================================================================
2622 LinearOrQuadratic::LinearOrQuadratic()
2627 void LinearOrQuadratic::SetMesh( const SMDS_Mesh* theMesh )
2632 bool LinearOrQuadratic::IsSatisfy( long theId )
2634 if (!myMesh) return false;
2635 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2636 if ( !anElem || (myType != SMDSAbs_All && anElem->GetType() != myType) )
2638 return (!anElem->IsQuadratic());
2641 void LinearOrQuadratic::SetType( SMDSAbs_ElementType theType )
2646 SMDSAbs_ElementType LinearOrQuadratic::GetType() const
2651 //================================================================================
2654 Description : Functor for check color of group to whic mesh element belongs to
2656 //================================================================================
2658 GroupColor::GroupColor()
2662 bool GroupColor::IsSatisfy( long theId )
2664 return (myIDs.find( theId ) != myIDs.end());
2667 void GroupColor::SetType( SMDSAbs_ElementType theType )
2672 SMDSAbs_ElementType GroupColor::GetType() const
2677 static bool isEqual( const Quantity_Color& theColor1,
2678 const Quantity_Color& theColor2 )
2680 // tolerance to compare colors
2681 const double tol = 5*1e-3;
2682 return ( fabs( theColor1.Red() - theColor2.Red() ) < tol &&
2683 fabs( theColor1.Green() - theColor2.Green() ) < tol &&
2684 fabs( theColor1.Blue() - theColor2.Blue() ) < tol );
2687 void GroupColor::SetMesh( const SMDS_Mesh* theMesh )
2691 const SMESHDS_Mesh* aMesh = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
2695 int nbGrp = aMesh->GetNbGroups();
2699 // iterates on groups and find necessary elements ids
2700 const std::set<SMESHDS_GroupBase*>& aGroups = aMesh->GetGroups();
2701 set<SMESHDS_GroupBase*>::const_iterator GrIt = aGroups.begin();
2702 for (; GrIt != aGroups.end(); GrIt++) {
2703 SMESHDS_GroupBase* aGrp = (*GrIt);
2706 // check type and color of group
2707 if ( !isEqual( myColor, aGrp->GetColor() ) )
2709 if ( myType != SMDSAbs_All && myType != (SMDSAbs_ElementType)aGrp->GetType() )
2712 SMDSAbs_ElementType aGrpElType = (SMDSAbs_ElementType)aGrp->GetType();
2713 if ( myType == aGrpElType || (myType == SMDSAbs_All && aGrpElType != SMDSAbs_Node) ) {
2714 // add elements IDS into control
2715 int aSize = aGrp->Extent();
2716 for (int i = 0; i < aSize; i++)
2717 myIDs.insert( aGrp->GetID(i+1) );
2722 void GroupColor::SetColorStr( const TCollection_AsciiString& theStr )
2724 Kernel_Utils::Localizer loc;
2725 TCollection_AsciiString aStr = theStr;
2726 aStr.RemoveAll( ' ' );
2727 aStr.RemoveAll( '\t' );
2728 for ( int aPos = aStr.Search( ";;" ); aPos != -1; aPos = aStr.Search( ";;" ) )
2729 aStr.Remove( aPos, 2 );
2730 Standard_Real clr[3];
2731 clr[0] = clr[1] = clr[2] = 0.;
2732 for ( int i = 0; i < 3; i++ ) {
2733 TCollection_AsciiString tmpStr = aStr.Token( ";", i+1 );
2734 if ( !tmpStr.IsEmpty() && tmpStr.IsRealValue() )
2735 clr[i] = tmpStr.RealValue();
2737 myColor = Quantity_Color( clr[0], clr[1], clr[2], Quantity_TOC_RGB );
2740 //=======================================================================
2741 // name : GetRangeStr
2742 // Purpose : Get range as a string.
2743 // Example: "1,2,3,50-60,63,67,70-"
2744 //=======================================================================
2746 void GroupColor::GetColorStr( TCollection_AsciiString& theResStr ) const
2749 theResStr += TCollection_AsciiString( myColor.Red() );
2750 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Green() );
2751 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Blue() );
2754 //================================================================================
2756 Class : ElemGeomType
2757 Description : Predicate to check element geometry type
2759 //================================================================================
2761 ElemGeomType::ElemGeomType()
2764 myType = SMDSAbs_All;
2765 myGeomType = SMDSGeom_TRIANGLE;
2768 void ElemGeomType::SetMesh( const SMDS_Mesh* theMesh )
2773 bool ElemGeomType::IsSatisfy( long theId )
2775 if (!myMesh) return false;
2776 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2779 const SMDSAbs_ElementType anElemType = anElem->GetType();
2780 if ( myType != SMDSAbs_All && anElemType != myType )
2782 bool isOk = ( anElem->GetGeomType() == myGeomType );
2786 void ElemGeomType::SetType( SMDSAbs_ElementType theType )
2791 SMDSAbs_ElementType ElemGeomType::GetType() const
2796 void ElemGeomType::SetGeomType( SMDSAbs_GeometryType theType )
2798 myGeomType = theType;
2801 SMDSAbs_GeometryType ElemGeomType::GetGeomType() const
2806 //================================================================================
2808 Class : ElemEntityType
2809 Description : Predicate to check element entity type
2811 //================================================================================
2813 ElemEntityType::ElemEntityType():
2815 myType( SMDSAbs_All ),
2816 myEntityType( SMDSEntity_0D )
2820 void ElemEntityType::SetMesh( const SMDS_Mesh* theMesh )
2825 bool ElemEntityType::IsSatisfy( long theId )
2827 if ( !myMesh ) return false;
2828 if ( myType == SMDSAbs_Node )
2829 return myMesh->FindNode( theId );
2830 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2832 myEntityType == anElem->GetEntityType() );
2835 void ElemEntityType::SetType( SMDSAbs_ElementType theType )
2840 SMDSAbs_ElementType ElemEntityType::GetType() const
2845 void ElemEntityType::SetElemEntityType( SMDSAbs_EntityType theEntityType )
2847 myEntityType = theEntityType;
2850 SMDSAbs_EntityType ElemEntityType::GetElemEntityType() const
2852 return myEntityType;
2855 //================================================================================
2857 * \brief Class ConnectedElements
2859 //================================================================================
2861 ConnectedElements::ConnectedElements():
2862 myNodeID(0), myType( SMDSAbs_All ), myOkIDsReady( false ) {}
2864 SMDSAbs_ElementType ConnectedElements::GetType() const
2867 int ConnectedElements::GetNode() const
2868 { return myXYZ.empty() ? myNodeID : 0; } // myNodeID can be found by myXYZ
2870 std::vector<double> ConnectedElements::GetPoint() const
2873 void ConnectedElements::clearOkIDs()
2874 { myOkIDsReady = false; myOkIDs.clear(); }
2876 void ConnectedElements::SetType( SMDSAbs_ElementType theType )
2878 if ( myType != theType || myMeshModifTracer.IsMeshModified() )
2883 void ConnectedElements::SetMesh( const SMDS_Mesh* theMesh )
2885 myMeshModifTracer.SetMesh( theMesh );
2886 if ( myMeshModifTracer.IsMeshModified() )
2889 if ( !myXYZ.empty() )
2890 SetPoint( myXYZ[0], myXYZ[1], myXYZ[2] ); // find a node near myXYZ it in a new mesh
2894 void ConnectedElements::SetNode( int nodeID )
2899 bool isSameDomain = false;
2900 if ( myOkIDsReady && myMeshModifTracer.GetMesh() && !myMeshModifTracer.IsMeshModified() )
2901 if ( const SMDS_MeshNode* n = myMeshModifTracer.GetMesh()->FindNode( myNodeID ))
2903 SMDS_ElemIteratorPtr eIt = n->GetInverseElementIterator( myType );
2904 while ( !isSameDomain && eIt->more() )
2905 isSameDomain = IsSatisfy( eIt->next()->GetID() );
2907 if ( !isSameDomain )
2911 void ConnectedElements::SetPoint( double x, double y, double z )
2919 bool isSameDomain = false;
2921 // find myNodeID by myXYZ if possible
2922 if ( myMeshModifTracer.GetMesh() )
2924 auto_ptr<SMESH_ElementSearcher> searcher
2925 ( SMESH_MeshAlgos::GetElementSearcher( (SMDS_Mesh&) *myMeshModifTracer.GetMesh() ));
2927 vector< const SMDS_MeshElement* > foundElems;
2928 searcher->FindElementsByPoint( gp_Pnt(x,y,z), SMDSAbs_All, foundElems );
2930 if ( !foundElems.empty() )
2932 myNodeID = foundElems[0]->GetNode(0)->GetID();
2933 if ( myOkIDsReady && !myMeshModifTracer.IsMeshModified() )
2934 isSameDomain = IsSatisfy( foundElems[0]->GetID() );
2937 if ( !isSameDomain )
2941 bool ConnectedElements::IsSatisfy( long theElementId )
2943 // Here we do NOT check if the mesh has changed, we do it in Set...() only!!!
2945 if ( !myOkIDsReady )
2947 if ( !myMeshModifTracer.GetMesh() )
2949 const SMDS_MeshNode* node0 = myMeshModifTracer.GetMesh()->FindNode( myNodeID );
2953 list< const SMDS_MeshNode* > nodeQueue( 1, node0 );
2954 std::set< int > checkedNodeIDs;
2956 // foreach node in nodeQueue:
2957 // foreach element sharing a node:
2958 // add ID of an element of myType to myOkIDs;
2959 // push all element nodes absent from checkedNodeIDs to nodeQueue;
2960 while ( !nodeQueue.empty() )
2962 const SMDS_MeshNode* node = nodeQueue.front();
2963 nodeQueue.pop_front();
2965 // loop on elements sharing the node
2966 SMDS_ElemIteratorPtr eIt = node->GetInverseElementIterator();
2967 while ( eIt->more() )
2969 // keep elements of myType
2970 const SMDS_MeshElement* element = eIt->next();
2971 if ( element->GetType() == myType )
2972 myOkIDs.insert( myOkIDs.end(), element->GetID() );
2974 // enqueue nodes of the element
2975 SMDS_ElemIteratorPtr nIt = element->nodesIterator();
2976 while ( nIt->more() )
2978 const SMDS_MeshNode* n = static_cast< const SMDS_MeshNode* >( nIt->next() );
2979 if ( checkedNodeIDs.insert( n->GetID() ).second )
2980 nodeQueue.push_back( n );
2984 if ( myType == SMDSAbs_Node )
2985 std::swap( myOkIDs, checkedNodeIDs );
2987 size_t totalNbElems = myMeshModifTracer.GetMesh()->GetMeshInfo().NbElements( myType );
2988 if ( myOkIDs.size() == totalNbElems )
2991 myOkIDsReady = true;
2994 return myOkIDs.empty() ? true : myOkIDs.count( theElementId );
2997 //================================================================================
2999 * \brief Class CoplanarFaces
3001 //================================================================================
3003 CoplanarFaces::CoplanarFaces()
3004 : myFaceID(0), myToler(0)
3007 void CoplanarFaces::SetMesh( const SMDS_Mesh* theMesh )
3009 myMeshModifTracer.SetMesh( theMesh );
3010 if ( myMeshModifTracer.IsMeshModified() )
3012 // Build a set of coplanar face ids
3014 myCoplanarIDs.clear();
3016 if ( !myMeshModifTracer.GetMesh() || !myFaceID || !myToler )
3019 const SMDS_MeshElement* face = myMeshModifTracer.GetMesh()->FindElement( myFaceID );
3020 if ( !face || face->GetType() != SMDSAbs_Face )
3024 gp_Vec myNorm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
3028 const double radianTol = myToler * M_PI / 180.;
3029 std::set< SMESH_TLink > checkedLinks;
3031 std::list< pair< const SMDS_MeshElement*, gp_Vec > > faceQueue;
3032 faceQueue.push_back( make_pair( face, myNorm ));
3033 while ( !faceQueue.empty() )
3035 face = faceQueue.front().first;
3036 myNorm = faceQueue.front().second;
3037 faceQueue.pop_front();
3039 for ( int i = 0, nbN = face->NbCornerNodes(); i < nbN; ++i )
3041 const SMDS_MeshNode* n1 = face->GetNode( i );
3042 const SMDS_MeshNode* n2 = face->GetNode(( i+1 )%nbN);
3043 if ( !checkedLinks.insert( SMESH_TLink( n1, n2 )).second )
3045 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator(SMDSAbs_Face);
3046 while ( fIt->more() )
3048 const SMDS_MeshElement* f = fIt->next();
3049 if ( f->GetNodeIndex( n2 ) > -1 )
3051 gp_Vec norm = getNormale( static_cast<const SMDS_MeshFace*>(f), &normOK );
3052 if (!normOK || myNorm.Angle( norm ) <= radianTol)
3054 myCoplanarIDs.insert( f->GetID() );
3055 faceQueue.push_back( make_pair( f, norm ));
3063 bool CoplanarFaces::IsSatisfy( long theElementId )
3065 return myCoplanarIDs.count( theElementId );
3070 *Description : Predicate for Range of Ids.
3071 * Range may be specified with two ways.
3072 * 1. Using AddToRange method
3073 * 2. With SetRangeStr method. Parameter of this method is a string
3074 * like as "1,2,3,50-60,63,67,70-"
3077 //=======================================================================
3078 // name : RangeOfIds
3079 // Purpose : Constructor
3080 //=======================================================================
3081 RangeOfIds::RangeOfIds()
3084 myType = SMDSAbs_All;
3087 //=======================================================================
3089 // Purpose : Set mesh
3090 //=======================================================================
3091 void RangeOfIds::SetMesh( const SMDS_Mesh* theMesh )
3096 //=======================================================================
3097 // name : AddToRange
3098 // Purpose : Add ID to the range
3099 //=======================================================================
3100 bool RangeOfIds::AddToRange( long theEntityId )
3102 myIds.Add( theEntityId );
3106 //=======================================================================
3107 // name : GetRangeStr
3108 // Purpose : Get range as a string.
3109 // Example: "1,2,3,50-60,63,67,70-"
3110 //=======================================================================
3111 void RangeOfIds::GetRangeStr( TCollection_AsciiString& theResStr )
3115 TColStd_SequenceOfInteger anIntSeq;
3116 TColStd_SequenceOfAsciiString aStrSeq;
3118 TColStd_MapIteratorOfMapOfInteger anIter( myIds );
3119 for ( ; anIter.More(); anIter.Next() )
3121 int anId = anIter.Key();
3122 TCollection_AsciiString aStr( anId );
3123 anIntSeq.Append( anId );
3124 aStrSeq.Append( aStr );
3127 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
3129 int aMinId = myMin( i );
3130 int aMaxId = myMax( i );
3132 TCollection_AsciiString aStr;
3133 if ( aMinId != IntegerFirst() )
3138 if ( aMaxId != IntegerLast() )
3141 // find position of the string in result sequence and insert string in it
3142 if ( anIntSeq.Length() == 0 )
3144 anIntSeq.Append( aMinId );
3145 aStrSeq.Append( aStr );
3149 if ( aMinId < anIntSeq.First() )
3151 anIntSeq.Prepend( aMinId );
3152 aStrSeq.Prepend( aStr );
3154 else if ( aMinId > anIntSeq.Last() )
3156 anIntSeq.Append( aMinId );
3157 aStrSeq.Append( aStr );
3160 for ( int j = 1, k = anIntSeq.Length(); j <= k; j++ )
3161 if ( aMinId < anIntSeq( j ) )
3163 anIntSeq.InsertBefore( j, aMinId );
3164 aStrSeq.InsertBefore( j, aStr );
3170 if ( aStrSeq.Length() == 0 )
3173 theResStr = aStrSeq( 1 );
3174 for ( int j = 2, k = aStrSeq.Length(); j <= k; j++ )
3177 theResStr += aStrSeq( j );
3181 //=======================================================================
3182 // name : SetRangeStr
3183 // Purpose : Define range with string
3184 // Example of entry string: "1,2,3,50-60,63,67,70-"
3185 //=======================================================================
3186 bool RangeOfIds::SetRangeStr( const TCollection_AsciiString& theStr )
3192 TCollection_AsciiString aStr = theStr;
3193 //aStr.RemoveAll( ' ' );
3194 //aStr.RemoveAll( '\t' );
3195 for ( int i = 1; i <= aStr.Length(); ++i )
3196 if ( isspace( aStr.Value( i )))
3197 aStr.SetValue( i, ',');
3199 for ( int aPos = aStr.Search( ",," ); aPos != -1; aPos = aStr.Search( ",," ) )
3200 aStr.Remove( aPos, 1 );
3202 TCollection_AsciiString tmpStr = aStr.Token( ",", 1 );
3204 while ( tmpStr != "" )
3206 tmpStr = aStr.Token( ",", i++ );
3207 int aPos = tmpStr.Search( '-' );
3211 if ( tmpStr.IsIntegerValue() )
3212 myIds.Add( tmpStr.IntegerValue() );
3218 TCollection_AsciiString aMaxStr = tmpStr.Split( aPos );
3219 TCollection_AsciiString aMinStr = tmpStr;
3221 while ( aMinStr.Search( "-" ) != -1 ) aMinStr.RemoveAll( '-' );
3222 while ( aMaxStr.Search( "-" ) != -1 ) aMaxStr.RemoveAll( '-' );
3224 if ( (!aMinStr.IsEmpty() && !aMinStr.IsIntegerValue()) ||
3225 (!aMaxStr.IsEmpty() && !aMaxStr.IsIntegerValue()) )
3228 myMin.Append( aMinStr.IsEmpty() ? IntegerFirst() : aMinStr.IntegerValue() );
3229 myMax.Append( aMaxStr.IsEmpty() ? IntegerLast() : aMaxStr.IntegerValue() );
3236 //=======================================================================
3238 // Purpose : Get type of supported entities
3239 //=======================================================================
3240 SMDSAbs_ElementType RangeOfIds::GetType() const
3245 //=======================================================================
3247 // Purpose : Set type of supported entities
3248 //=======================================================================
3249 void RangeOfIds::SetType( SMDSAbs_ElementType theType )
3254 //=======================================================================
3256 // Purpose : Verify whether entity satisfies to this rpedicate
3257 //=======================================================================
3258 bool RangeOfIds::IsSatisfy( long theId )
3263 if ( myType == SMDSAbs_Node )
3265 if ( myMesh->FindNode( theId ) == 0 )
3270 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
3271 if ( anElem == 0 || (myType != anElem->GetType() && myType != SMDSAbs_All ))
3275 if ( myIds.Contains( theId ) )
3278 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
3279 if ( theId >= myMin( i ) && theId <= myMax( i ) )
3287 Description : Base class for comparators
3289 Comparator::Comparator():
3293 Comparator::~Comparator()
3296 void Comparator::SetMesh( const SMDS_Mesh* theMesh )
3299 myFunctor->SetMesh( theMesh );
3302 void Comparator::SetMargin( double theValue )
3304 myMargin = theValue;
3307 void Comparator::SetNumFunctor( NumericalFunctorPtr theFunct )
3309 myFunctor = theFunct;
3312 SMDSAbs_ElementType Comparator::GetType() const
3314 return myFunctor ? myFunctor->GetType() : SMDSAbs_All;
3317 double Comparator::GetMargin()
3325 Description : Comparator "<"
3327 bool LessThan::IsSatisfy( long theId )
3329 return myFunctor && myFunctor->GetValue( theId ) < myMargin;
3335 Description : Comparator ">"
3337 bool MoreThan::IsSatisfy( long theId )
3339 return myFunctor && myFunctor->GetValue( theId ) > myMargin;
3345 Description : Comparator "="
3348 myToler(Precision::Confusion())
3351 bool EqualTo::IsSatisfy( long theId )
3353 return myFunctor && fabs( myFunctor->GetValue( theId ) - myMargin ) < myToler;
3356 void EqualTo::SetTolerance( double theToler )
3361 double EqualTo::GetTolerance()
3368 Description : Logical NOT predicate
3370 LogicalNOT::LogicalNOT()
3373 LogicalNOT::~LogicalNOT()
3376 bool LogicalNOT::IsSatisfy( long theId )
3378 return myPredicate && !myPredicate->IsSatisfy( theId );
3381 void LogicalNOT::SetMesh( const SMDS_Mesh* theMesh )
3384 myPredicate->SetMesh( theMesh );
3387 void LogicalNOT::SetPredicate( PredicatePtr thePred )
3389 myPredicate = thePred;
3392 SMDSAbs_ElementType LogicalNOT::GetType() const
3394 return myPredicate ? myPredicate->GetType() : SMDSAbs_All;
3399 Class : LogicalBinary
3400 Description : Base class for binary logical predicate
3402 LogicalBinary::LogicalBinary()
3405 LogicalBinary::~LogicalBinary()
3408 void LogicalBinary::SetMesh( const SMDS_Mesh* theMesh )
3411 myPredicate1->SetMesh( theMesh );
3414 myPredicate2->SetMesh( theMesh );
3417 void LogicalBinary::SetPredicate1( PredicatePtr thePredicate )
3419 myPredicate1 = thePredicate;
3422 void LogicalBinary::SetPredicate2( PredicatePtr thePredicate )
3424 myPredicate2 = thePredicate;
3427 SMDSAbs_ElementType LogicalBinary::GetType() const
3429 if ( !myPredicate1 || !myPredicate2 )
3432 SMDSAbs_ElementType aType1 = myPredicate1->GetType();
3433 SMDSAbs_ElementType aType2 = myPredicate2->GetType();
3435 return aType1 == aType2 ? aType1 : SMDSAbs_All;
3441 Description : Logical AND
3443 bool LogicalAND::IsSatisfy( long theId )
3448 myPredicate1->IsSatisfy( theId ) &&
3449 myPredicate2->IsSatisfy( theId );
3455 Description : Logical OR
3457 bool LogicalOR::IsSatisfy( long theId )
3462 (myPredicate1->IsSatisfy( theId ) ||
3463 myPredicate2->IsSatisfy( theId ));
3472 // #include <tbb/parallel_for.h>
3473 // #include <tbb/enumerable_thread_specific.h>
3475 // namespace Parallel
3477 // typedef tbb::enumerable_thread_specific< TIdSequence > TIdSeq;
3481 // const SMDS_Mesh* myMesh;
3482 // PredicatePtr myPredicate;
3483 // TIdSeq & myOKIds;
3484 // Predicate( const SMDS_Mesh* m, PredicatePtr p, TIdSeq & ids ):
3485 // myMesh(m), myPredicate(p->Duplicate()), myOKIds(ids) {}
3486 // void operator() ( const tbb::blocked_range<size_t>& r ) const
3488 // for ( size_t i = r.begin(); i != r.end(); ++i )
3489 // if ( myPredicate->IsSatisfy( i ))
3490 // myOKIds.local().push_back();
3502 void Filter::SetPredicate( PredicatePtr thePredicate )
3504 myPredicate = thePredicate;
3507 void Filter::GetElementsId( const SMDS_Mesh* theMesh,
3508 PredicatePtr thePredicate,
3509 TIdSequence& theSequence )
3511 theSequence.clear();
3513 if ( !theMesh || !thePredicate )
3516 thePredicate->SetMesh( theMesh );
3518 SMDS_ElemIteratorPtr elemIt = theMesh->elementsIterator( thePredicate->GetType() );
3520 while ( elemIt->more() ) {
3521 const SMDS_MeshElement* anElem = elemIt->next();
3522 long anId = anElem->GetID();
3523 if ( thePredicate->IsSatisfy( anId ) )
3524 theSequence.push_back( anId );
3529 void Filter::GetElementsId( const SMDS_Mesh* theMesh,
3530 Filter::TIdSequence& theSequence )
3532 GetElementsId(theMesh,myPredicate,theSequence);
3539 typedef std::set<SMDS_MeshFace*> TMapOfFacePtr;
3545 ManifoldPart::Link::Link( SMDS_MeshNode* theNode1,
3546 SMDS_MeshNode* theNode2 )
3552 ManifoldPart::Link::~Link()
3558 bool ManifoldPart::Link::IsEqual( const ManifoldPart::Link& theLink ) const
3560 if ( myNode1 == theLink.myNode1 &&
3561 myNode2 == theLink.myNode2 )
3563 else if ( myNode1 == theLink.myNode2 &&
3564 myNode2 == theLink.myNode1 )
3570 bool ManifoldPart::Link::operator<( const ManifoldPart::Link& x ) const
3572 if(myNode1 < x.myNode1) return true;
3573 if(myNode1 == x.myNode1)
3574 if(myNode2 < x.myNode2) return true;
3578 bool ManifoldPart::IsEqual( const ManifoldPart::Link& theLink1,
3579 const ManifoldPart::Link& theLink2 )
3581 return theLink1.IsEqual( theLink2 );
3584 ManifoldPart::ManifoldPart()
3587 myAngToler = Precision::Angular();
3588 myIsOnlyManifold = true;
3591 ManifoldPart::~ManifoldPart()
3596 void ManifoldPart::SetMesh( const SMDS_Mesh* theMesh )
3602 SMDSAbs_ElementType ManifoldPart::GetType() const
3603 { return SMDSAbs_Face; }
3605 bool ManifoldPart::IsSatisfy( long theElementId )
3607 return myMapIds.Contains( theElementId );
3610 void ManifoldPart::SetAngleTolerance( const double theAngToler )
3611 { myAngToler = theAngToler; }
3613 double ManifoldPart::GetAngleTolerance() const
3614 { return myAngToler; }
3616 void ManifoldPart::SetIsOnlyManifold( const bool theIsOnly )
3617 { myIsOnlyManifold = theIsOnly; }
3619 void ManifoldPart::SetStartElem( const long theStartId )
3620 { myStartElemId = theStartId; }
3622 bool ManifoldPart::process()
3625 myMapBadGeomIds.Clear();
3627 myAllFacePtr.clear();
3628 myAllFacePtrIntDMap.clear();
3632 // collect all faces into own map
3633 SMDS_FaceIteratorPtr anFaceItr = myMesh->facesIterator();
3634 for (; anFaceItr->more(); )
3636 SMDS_MeshFace* aFacePtr = (SMDS_MeshFace*)anFaceItr->next();
3637 myAllFacePtr.push_back( aFacePtr );
3638 myAllFacePtrIntDMap[aFacePtr] = myAllFacePtr.size()-1;
3641 SMDS_MeshFace* aStartFace = (SMDS_MeshFace*)myMesh->FindElement( myStartElemId );
3645 // the map of non manifold links and bad geometry
3646 TMapOfLink aMapOfNonManifold;
3647 TColStd_MapOfInteger aMapOfTreated;
3649 // begin cycle on faces from start index and run on vector till the end
3650 // and from begin to start index to cover whole vector
3651 const int aStartIndx = myAllFacePtrIntDMap[aStartFace];
3652 bool isStartTreat = false;
3653 for ( int fi = aStartIndx; !isStartTreat || fi != aStartIndx ; fi++ )
3655 if ( fi == aStartIndx )
3656 isStartTreat = true;
3657 // as result next time when fi will be equal to aStartIndx
3659 SMDS_MeshFace* aFacePtr = myAllFacePtr[ fi ];
3660 if ( aMapOfTreated.Contains( aFacePtr->GetID() ) )
3663 aMapOfTreated.Add( aFacePtr->GetID() );
3664 TColStd_MapOfInteger aResFaces;
3665 if ( !findConnected( myAllFacePtrIntDMap, aFacePtr,
3666 aMapOfNonManifold, aResFaces ) )
3668 TColStd_MapIteratorOfMapOfInteger anItr( aResFaces );
3669 for ( ; anItr.More(); anItr.Next() )
3671 int aFaceId = anItr.Key();
3672 aMapOfTreated.Add( aFaceId );
3673 myMapIds.Add( aFaceId );
3676 if ( fi == ( myAllFacePtr.size() - 1 ) )
3678 } // end run on vector of faces
3679 return !myMapIds.IsEmpty();
3682 static void getLinks( const SMDS_MeshFace* theFace,
3683 ManifoldPart::TVectorOfLink& theLinks )
3685 int aNbNode = theFace->NbNodes();
3686 SMDS_ElemIteratorPtr aNodeItr = theFace->nodesIterator();
3688 SMDS_MeshNode* aNode = 0;
3689 for ( ; aNodeItr->more() && i <= aNbNode; )
3692 SMDS_MeshNode* aN1 = (SMDS_MeshNode*)aNodeItr->next();
3696 SMDS_MeshNode* aN2 = ( i >= aNbNode ) ? aNode : (SMDS_MeshNode*)aNodeItr->next();
3698 ManifoldPart::Link aLink( aN1, aN2 );
3699 theLinks.push_back( aLink );
3703 bool ManifoldPart::findConnected
3704 ( const ManifoldPart::TDataMapFacePtrInt& theAllFacePtrInt,
3705 SMDS_MeshFace* theStartFace,
3706 ManifoldPart::TMapOfLink& theNonManifold,
3707 TColStd_MapOfInteger& theResFaces )
3709 theResFaces.Clear();
3710 if ( !theAllFacePtrInt.size() )
3713 if ( getNormale( theStartFace ).SquareModulus() <= gp::Resolution() )
3715 myMapBadGeomIds.Add( theStartFace->GetID() );
3719 ManifoldPart::TMapOfLink aMapOfBoundary, aMapToSkip;
3720 ManifoldPart::TVectorOfLink aSeqOfBoundary;
3721 theResFaces.Add( theStartFace->GetID() );
3722 ManifoldPart::TDataMapOfLinkFacePtr aDMapLinkFace;
3724 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3725 aDMapLinkFace, theNonManifold, theStartFace );
3727 bool isDone = false;
3728 while ( !isDone && aMapOfBoundary.size() != 0 )
3730 bool isToReset = false;
3731 ManifoldPart::TVectorOfLink::iterator pLink = aSeqOfBoundary.begin();
3732 for ( ; !isToReset && pLink != aSeqOfBoundary.end(); ++pLink )
3734 ManifoldPart::Link aLink = *pLink;
3735 if ( aMapToSkip.find( aLink ) != aMapToSkip.end() )
3737 // each link could be treated only once
3738 aMapToSkip.insert( aLink );
3740 ManifoldPart::TVectorOfFacePtr aFaces;
3742 if ( myIsOnlyManifold &&
3743 (theNonManifold.find( aLink ) != theNonManifold.end()) )
3747 getFacesByLink( aLink, aFaces );
3748 // filter the element to keep only indicated elements
3749 ManifoldPart::TVectorOfFacePtr aFiltered;
3750 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3751 for ( ; pFace != aFaces.end(); ++pFace )
3753 SMDS_MeshFace* aFace = *pFace;
3754 if ( myAllFacePtrIntDMap.find( aFace ) != myAllFacePtrIntDMap.end() )
3755 aFiltered.push_back( aFace );
3758 if ( aFaces.size() < 2 ) // no neihgbour faces
3760 else if ( myIsOnlyManifold && aFaces.size() > 2 ) // non manifold case
3762 theNonManifold.insert( aLink );
3767 // compare normal with normals of neighbor element
3768 SMDS_MeshFace* aPrevFace = aDMapLinkFace[ aLink ];
3769 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3770 for ( ; pFace != aFaces.end(); ++pFace )
3772 SMDS_MeshFace* aNextFace = *pFace;
3773 if ( aPrevFace == aNextFace )
3775 int anNextFaceID = aNextFace->GetID();
3776 if ( myIsOnlyManifold && theResFaces.Contains( anNextFaceID ) )
3777 // should not be with non manifold restriction. probably bad topology
3779 // check if face was treated and skipped
3780 if ( myMapBadGeomIds.Contains( anNextFaceID ) ||
3781 !isInPlane( aPrevFace, aNextFace ) )
3783 // add new element to connected and extend the boundaries.
3784 theResFaces.Add( anNextFaceID );
3785 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3786 aDMapLinkFace, theNonManifold, aNextFace );
3790 isDone = !isToReset;
3793 return !theResFaces.IsEmpty();
3796 bool ManifoldPart::isInPlane( const SMDS_MeshFace* theFace1,
3797 const SMDS_MeshFace* theFace2 )
3799 gp_Dir aNorm1 = gp_Dir( getNormale( theFace1 ) );
3800 gp_XYZ aNorm2XYZ = getNormale( theFace2 );
3801 if ( aNorm2XYZ.SquareModulus() <= gp::Resolution() )
3803 myMapBadGeomIds.Add( theFace2->GetID() );
3806 if ( aNorm1.IsParallel( gp_Dir( aNorm2XYZ ), myAngToler ) )
3812 void ManifoldPart::expandBoundary
3813 ( ManifoldPart::TMapOfLink& theMapOfBoundary,
3814 ManifoldPart::TVectorOfLink& theSeqOfBoundary,
3815 ManifoldPart::TDataMapOfLinkFacePtr& theDMapLinkFacePtr,
3816 ManifoldPart::TMapOfLink& theNonManifold,
3817 SMDS_MeshFace* theNextFace ) const
3819 ManifoldPart::TVectorOfLink aLinks;
3820 getLinks( theNextFace, aLinks );
3821 int aNbLink = (int)aLinks.size();
3822 for ( int i = 0; i < aNbLink; i++ )
3824 ManifoldPart::Link aLink = aLinks[ i ];
3825 if ( myIsOnlyManifold && (theNonManifold.find( aLink ) != theNonManifold.end()) )
3827 if ( theMapOfBoundary.find( aLink ) != theMapOfBoundary.end() )
3829 if ( myIsOnlyManifold )
3831 // remove from boundary
3832 theMapOfBoundary.erase( aLink );
3833 ManifoldPart::TVectorOfLink::iterator pLink = theSeqOfBoundary.begin();
3834 for ( ; pLink != theSeqOfBoundary.end(); ++pLink )
3836 ManifoldPart::Link aBoundLink = *pLink;
3837 if ( aBoundLink.IsEqual( aLink ) )
3839 theSeqOfBoundary.erase( pLink );
3847 theMapOfBoundary.insert( aLink );
3848 theSeqOfBoundary.push_back( aLink );
3849 theDMapLinkFacePtr[ aLink ] = theNextFace;
3854 void ManifoldPart::getFacesByLink( const ManifoldPart::Link& theLink,
3855 ManifoldPart::TVectorOfFacePtr& theFaces ) const
3857 std::set<SMDS_MeshCell *> aSetOfFaces;
3858 // take all faces that shared first node
3859 SMDS_ElemIteratorPtr anItr = theLink.myNode1->facesIterator();
3860 for ( ; anItr->more(); )
3862 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3865 aSetOfFaces.insert( aFace );
3867 // take all faces that shared second node
3868 anItr = theLink.myNode2->facesIterator();
3869 // find the common part of two sets
3870 for ( ; anItr->more(); )
3872 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3873 if ( aSetOfFaces.count( aFace ) )
3874 theFaces.push_back( aFace );
3879 Class : BelongToMeshGroup
3880 Description : Verify whether a mesh element is included into a mesh group
3882 BelongToMeshGroup::BelongToMeshGroup(): myGroup( 0 )
3886 void BelongToMeshGroup::SetGroup( SMESHDS_GroupBase* g )
3891 void BelongToMeshGroup::SetStoreName( const std::string& sn )
3896 void BelongToMeshGroup::SetMesh( const SMDS_Mesh* theMesh )
3898 if ( myGroup && myGroup->GetMesh() != theMesh )
3902 if ( !myGroup && !myStoreName.empty() )
3904 if ( const SMESHDS_Mesh* aMesh = dynamic_cast<const SMESHDS_Mesh*>(theMesh))
3906 const std::set<SMESHDS_GroupBase*>& grps = aMesh->GetGroups();
3907 std::set<SMESHDS_GroupBase*>::const_iterator g = grps.begin();
3908 for ( ; g != grps.end() && !myGroup; ++g )
3909 if ( *g && myStoreName == (*g)->GetStoreName() )
3915 myGroup->IsEmpty(); // make GroupOnFilter update its predicate
3919 bool BelongToMeshGroup::IsSatisfy( long theElementId )
3921 return myGroup ? myGroup->Contains( theElementId ) : false;
3924 SMDSAbs_ElementType BelongToMeshGroup::GetType() const
3926 return myGroup ? myGroup->GetType() : SMDSAbs_All;
3933 ElementsOnSurface::ElementsOnSurface()
3936 myType = SMDSAbs_All;
3938 myToler = Precision::Confusion();
3939 myUseBoundaries = false;
3942 ElementsOnSurface::~ElementsOnSurface()
3946 void ElementsOnSurface::SetMesh( const SMDS_Mesh* theMesh )
3948 myMeshModifTracer.SetMesh( theMesh );
3949 if ( myMeshModifTracer.IsMeshModified())
3953 bool ElementsOnSurface::IsSatisfy( long theElementId )
3955 return myIds.Contains( theElementId );
3958 SMDSAbs_ElementType ElementsOnSurface::GetType() const
3961 void ElementsOnSurface::SetTolerance( const double theToler )
3963 if ( myToler != theToler )
3968 double ElementsOnSurface::GetTolerance() const
3971 void ElementsOnSurface::SetUseBoundaries( bool theUse )
3973 if ( myUseBoundaries != theUse ) {
3974 myUseBoundaries = theUse;
3975 SetSurface( mySurf, myType );
3979 void ElementsOnSurface::SetSurface( const TopoDS_Shape& theShape,
3980 const SMDSAbs_ElementType theType )
3985 if ( theShape.IsNull() || theShape.ShapeType() != TopAbs_FACE )
3987 mySurf = TopoDS::Face( theShape );
3988 BRepAdaptor_Surface SA( mySurf, myUseBoundaries );
3990 u1 = SA.FirstUParameter(),
3991 u2 = SA.LastUParameter(),
3992 v1 = SA.FirstVParameter(),
3993 v2 = SA.LastVParameter();
3994 Handle(Geom_Surface) surf = BRep_Tool::Surface( mySurf );
3995 myProjector.Init( surf, u1,u2, v1,v2 );
3999 void ElementsOnSurface::process()
4002 if ( mySurf.IsNull() )
4005 if ( !myMeshModifTracer.GetMesh() )
4008 myIds.ReSize( myMeshModifTracer.GetMesh()->GetMeshInfo().NbElements( myType ));
4010 SMDS_ElemIteratorPtr anIter = myMeshModifTracer.GetMesh()->elementsIterator( myType );
4011 for(; anIter->more(); )
4012 process( anIter->next() );
4015 void ElementsOnSurface::process( const SMDS_MeshElement* theElemPtr )
4017 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
4018 bool isSatisfy = true;
4019 for ( ; aNodeItr->more(); )
4021 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
4022 if ( !isOnSurface( aNode ) )
4029 myIds.Add( theElemPtr->GetID() );
4032 bool ElementsOnSurface::isOnSurface( const SMDS_MeshNode* theNode )
4034 if ( mySurf.IsNull() )
4037 gp_Pnt aPnt( theNode->X(), theNode->Y(), theNode->Z() );
4038 // double aToler2 = myToler * myToler;
4039 // if ( mySurf->IsKind(STANDARD_TYPE(Geom_Plane)))
4041 // gp_Pln aPln = Handle(Geom_Plane)::DownCast(mySurf)->Pln();
4042 // if ( aPln.SquareDistance( aPnt ) > aToler2 )
4045 // else if ( mySurf->IsKind(STANDARD_TYPE(Geom_CylindricalSurface)))
4047 // gp_Cylinder aCyl = Handle(Geom_CylindricalSurface)::DownCast(mySurf)->Cylinder();
4048 // double aRad = aCyl.Radius();
4049 // gp_Ax3 anAxis = aCyl.Position();
4050 // gp_XYZ aLoc = aCyl.Location().XYZ();
4051 // double aXDist = anAxis.XDirection().XYZ() * ( aPnt.XYZ() - aLoc );
4052 // double aYDist = anAxis.YDirection().XYZ() * ( aPnt.XYZ() - aLoc );
4053 // if ( fabs(aXDist*aXDist + aYDist*aYDist - aRad*aRad) > aToler2 )
4058 myProjector.Perform( aPnt );
4059 bool isOn = ( myProjector.IsDone() && myProjector.LowerDistance() <= myToler );
4069 ElementsOnShape::ElementsOnShape()
4071 myType(SMDSAbs_All),
4072 myToler(Precision::Confusion()),
4073 myAllNodesFlag(false)
4077 ElementsOnShape::~ElementsOnShape()
4082 SMDSAbs_ElementType ElementsOnShape::GetType() const
4087 void ElementsOnShape::SetTolerance (const double theToler)
4089 if (myToler != theToler) {
4091 SetShape(myShape, myType);
4095 double ElementsOnShape::GetTolerance() const
4100 void ElementsOnShape::SetAllNodes (bool theAllNodes)
4102 myAllNodesFlag = theAllNodes;
4105 void ElementsOnShape::SetMesh (const SMDS_Mesh* theMesh)
4107 myMeshModifTracer.SetMesh( theMesh );
4108 if ( myMeshModifTracer.IsMeshModified())
4110 size_t nbNodes = theMesh ? theMesh->NbNodes() : 0;
4111 if ( myNodeIsChecked.size() == nbNodes )
4113 std::fill( myNodeIsChecked.begin(), myNodeIsChecked.end(), false );
4117 SMESHUtils::FreeVector( myNodeIsChecked );
4118 SMESHUtils::FreeVector( myNodeIsOut );
4119 myNodeIsChecked.resize( nbNodes, false );
4120 myNodeIsOut.resize( nbNodes );
4125 bool ElementsOnShape::getNodeIsOut( const SMDS_MeshNode* n, bool& isOut )
4127 if ( n->GetID() >= (int) myNodeIsChecked.size() ||
4128 !myNodeIsChecked[ n->GetID() ])
4131 isOut = myNodeIsOut[ n->GetID() ];
4135 void ElementsOnShape::setNodeIsOut( const SMDS_MeshNode* n, bool isOut )
4137 if ( n->GetID() < (int) myNodeIsChecked.size() )
4139 myNodeIsChecked[ n->GetID() ] = true;
4140 myNodeIsOut [ n->GetID() ] = isOut;
4144 void ElementsOnShape::SetShape (const TopoDS_Shape& theShape,
4145 const SMDSAbs_ElementType theType)
4149 if ( myShape.IsNull() ) return;
4151 TopTools_IndexedMapOfShape shapesMap;
4152 TopAbs_ShapeEnum shapeTypes[4] = { TopAbs_SOLID, TopAbs_FACE, TopAbs_EDGE, TopAbs_VERTEX };
4153 TopExp_Explorer sub;
4154 for ( int i = 0; i < 4; ++i )
4156 if ( shapesMap.IsEmpty() )
4157 for ( sub.Init( myShape, shapeTypes[i] ); sub.More(); sub.Next() )
4158 shapesMap.Add( sub.Current() );
4160 for ( sub.Init( myShape, shapeTypes[i], shapeTypes[i-1] ); sub.More(); sub.Next() )
4161 shapesMap.Add( sub.Current() );
4165 myClassifiers.resize( shapesMap.Extent() );
4166 for ( int i = 0; i < shapesMap.Extent(); ++i )
4167 myClassifiers[ i ] = new TClassifier( shapesMap( i+1 ), myToler );
4169 if ( theType == SMDSAbs_Node )
4171 SMESHUtils::FreeVector( myNodeIsChecked );
4172 SMESHUtils::FreeVector( myNodeIsOut );
4176 std::fill( myNodeIsChecked.begin(), myNodeIsChecked.end(), false );
4180 void ElementsOnShape::clearClassifiers()
4182 for ( size_t i = 0; i < myClassifiers.size(); ++i )
4183 delete myClassifiers[ i ];
4184 myClassifiers.clear();
4187 bool ElementsOnShape::IsSatisfy (long elemId)
4189 const SMDS_Mesh* mesh = myMeshModifTracer.GetMesh();
4190 const SMDS_MeshElement* elem =
4191 ( myType == SMDSAbs_Node ? mesh->FindNode( elemId ) : mesh->FindElement( elemId ));
4192 if ( !elem || myClassifiers.empty() )
4195 bool isSatisfy = myAllNodesFlag, isNodeOut;
4197 gp_XYZ centerXYZ (0, 0, 0);
4199 SMDS_ElemIteratorPtr aNodeItr = elem->nodesIterator();
4200 while (aNodeItr->more() && (isSatisfy == myAllNodesFlag))
4202 SMESH_TNodeXYZ aPnt( aNodeItr->next() );
4206 if ( !getNodeIsOut( aPnt._node, isNodeOut ))
4208 for ( size_t i = 0; i < myClassifiers.size() && isNodeOut; ++i )
4209 isNodeOut = myClassifiers[i]->IsOut( aPnt );
4211 setNodeIsOut( aPnt._node, isNodeOut );
4213 isSatisfy = !isNodeOut;
4216 // Check the center point for volumes MantisBug 0020168
4219 myClassifiers[0]->ShapeType() == TopAbs_SOLID)
4221 centerXYZ /= elem->NbNodes();
4223 for ( size_t i = 0; i < myClassifiers.size() && !isSatisfy; ++i )
4224 isSatisfy = ! myClassifiers[i]->IsOut( centerXYZ );
4230 TopAbs_ShapeEnum ElementsOnShape::TClassifier::ShapeType() const
4232 return myShape.ShapeType();
4235 bool ElementsOnShape::TClassifier::IsOut(const gp_Pnt& p)
4237 return (this->*myIsOutFun)( p );
4240 void ElementsOnShape::TClassifier::Init (const TopoDS_Shape& theShape, double theTol)
4244 switch ( myShape.ShapeType() )
4246 case TopAbs_SOLID: {
4247 if ( isBox( theShape ))
4249 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfBox;
4253 mySolidClfr.Load(theShape);
4254 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfSolid;
4259 Standard_Real u1,u2,v1,v2;
4260 Handle(Geom_Surface) surf = BRep_Tool::Surface( TopoDS::Face( theShape ));
4261 surf->Bounds( u1,u2,v1,v2 );
4262 myProjFace.Init(surf, u1,u2, v1,v2, myTol );
4263 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfFace;
4267 Standard_Real u1, u2;
4268 Handle(Geom_Curve) curve = BRep_Tool::Curve( TopoDS::Edge(theShape), u1, u2);
4269 myProjEdge.Init(curve, u1, u2);
4270 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfEdge;
4273 case TopAbs_VERTEX:{
4274 myVertexXYZ = BRep_Tool::Pnt( TopoDS::Vertex( theShape ) );
4275 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfVertex;
4279 throw SALOME_Exception("Programmer error in usage of ElementsOnShape::TClassifier");
4283 bool ElementsOnShape::TClassifier::isOutOfSolid (const gp_Pnt& p)
4285 mySolidClfr.Perform( p, myTol );
4286 return ( mySolidClfr.State() != TopAbs_IN && mySolidClfr.State() != TopAbs_ON );
4289 bool ElementsOnShape::TClassifier::isOutOfBox (const gp_Pnt& p)
4291 return myBox.IsOut( p.XYZ() );
4294 bool ElementsOnShape::TClassifier::isOutOfFace (const gp_Pnt& p)
4296 myProjFace.Perform( p );
4297 if ( myProjFace.IsDone() && myProjFace.LowerDistance() <= myTol )
4299 // check relatively to the face
4300 Quantity_Parameter u, v;
4301 myProjFace.LowerDistanceParameters(u, v);
4302 gp_Pnt2d aProjPnt (u, v);
4303 BRepClass_FaceClassifier aClsf ( TopoDS::Face( myShape ), aProjPnt, myTol );
4304 if ( aClsf.State() == TopAbs_IN || aClsf.State() == TopAbs_ON )
4310 bool ElementsOnShape::TClassifier::isOutOfEdge (const gp_Pnt& p)
4312 myProjEdge.Perform( p );
4313 return ! ( myProjEdge.NbPoints() > 0 && myProjEdge.LowerDistance() <= myTol );
4316 bool ElementsOnShape::TClassifier::isOutOfVertex(const gp_Pnt& p)
4318 return ( myVertexXYZ.Distance( p ) > myTol );
4321 bool ElementsOnShape::TClassifier::isBox (const TopoDS_Shape& theShape)
4323 TopTools_IndexedMapOfShape vMap;
4324 TopExp::MapShapes( theShape, TopAbs_VERTEX, vMap );
4325 if ( vMap.Extent() != 8 )
4329 for ( int i = 1; i <= 8; ++i )
4330 myBox.Add( BRep_Tool::Pnt( TopoDS::Vertex( vMap( i ))).XYZ() );
4332 gp_XYZ pMin = myBox.CornerMin(), pMax = myBox.CornerMax();
4333 for ( int i = 1; i <= 8; ++i )
4335 gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex( vMap( i )));
4336 for ( int iC = 1; iC <= 3; ++ iC )
4338 double d1 = Abs( pMin.Coord( iC ) - p.Coord( iC ));
4339 double d2 = Abs( pMax.Coord( iC ) - p.Coord( iC ));
4340 if ( Min( d1, d2 ) > myTol )
4344 myBox.Enlarge( myTol );
4350 Class : BelongToGeom
4351 Description : Predicate for verifying whether entity belongs to
4352 specified geometrical support
4355 BelongToGeom::BelongToGeom()
4357 myType(SMDSAbs_All),
4358 myIsSubshape(false),
4359 myTolerance(Precision::Confusion())
4362 void BelongToGeom::SetMesh( const SMDS_Mesh* theMesh )
4364 myMeshDS = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
4368 void BelongToGeom::SetGeom( const TopoDS_Shape& theShape )
4374 static bool IsSubShape (const TopTools_IndexedMapOfShape& theMap,
4375 const TopoDS_Shape& theShape)
4377 if (theMap.Contains(theShape)) return true;
4379 if (theShape.ShapeType() == TopAbs_COMPOUND ||
4380 theShape.ShapeType() == TopAbs_COMPSOLID)
4382 TopoDS_Iterator anIt (theShape, Standard_True, Standard_True);
4383 for (; anIt.More(); anIt.Next())
4385 if (!IsSubShape(theMap, anIt.Value())) {
4395 void BelongToGeom::init()
4397 if (!myMeshDS || myShape.IsNull()) return;
4399 // is sub-shape of main shape?
4400 TopoDS_Shape aMainShape = myMeshDS->ShapeToMesh();
4401 if (aMainShape.IsNull()) {
4402 myIsSubshape = false;
4405 TopTools_IndexedMapOfShape aMap;
4406 TopExp::MapShapes(aMainShape, aMap);
4407 myIsSubshape = IsSubShape(aMap, myShape);
4410 //if (!myIsSubshape) // to be always ready to check an element not bound to geometry
4412 myElementsOnShapePtr.reset(new ElementsOnShape());
4413 myElementsOnShapePtr->SetTolerance(myTolerance);
4414 myElementsOnShapePtr->SetAllNodes(true); // "belong", while false means "lays on"
4415 myElementsOnShapePtr->SetMesh(myMeshDS);
4416 myElementsOnShapePtr->SetShape(myShape, myType);
4420 static bool IsContains( const SMESHDS_Mesh* theMeshDS,
4421 const TopoDS_Shape& theShape,
4422 const SMDS_MeshElement* theElem,
4423 TopAbs_ShapeEnum theFindShapeEnum,
4424 TopAbs_ShapeEnum theAvoidShapeEnum = TopAbs_SHAPE )
4426 TopExp_Explorer anExp( theShape,theFindShapeEnum,theAvoidShapeEnum );
4428 while( anExp.More() )
4430 const TopoDS_Shape& aShape = anExp.Current();
4431 if( SMESHDS_SubMesh* aSubMesh = theMeshDS->MeshElements( aShape ) ){
4432 if( aSubMesh->Contains( theElem ) )
4440 bool BelongToGeom::IsSatisfy (long theId)
4442 if (myMeshDS == 0 || myShape.IsNull())
4447 return myElementsOnShapePtr->IsSatisfy(theId);
4451 if (myType == SMDSAbs_Node)
4453 if( const SMDS_MeshNode* aNode = myMeshDS->FindNode( theId ) )
4455 if ( aNode->getshapeId() < 1 )
4456 return myElementsOnShapePtr->IsSatisfy(theId);
4458 const SMDS_PositionPtr& aPosition = aNode->GetPosition();
4459 SMDS_TypeOfPosition aTypeOfPosition = aPosition->GetTypeOfPosition();
4460 switch( aTypeOfPosition )
4462 case SMDS_TOP_VERTEX : return ( IsContains( myMeshDS,myShape,aNode,TopAbs_VERTEX ));
4463 case SMDS_TOP_EDGE : return ( IsContains( myMeshDS,myShape,aNode,TopAbs_EDGE ));
4464 case SMDS_TOP_FACE : return ( IsContains( myMeshDS,myShape,aNode,TopAbs_FACE ));
4465 case SMDS_TOP_3DSPACE: return ( IsContains( myMeshDS,myShape,aNode,TopAbs_SOLID ) ||
4466 IsContains( myMeshDS,myShape,aNode,TopAbs_SHELL ));
4472 if ( const SMDS_MeshElement* anElem = myMeshDS->FindElement( theId ))
4474 if ( anElem->getshapeId() < 1 )
4475 return myElementsOnShapePtr->IsSatisfy(theId);
4477 if( myType == SMDSAbs_All )
4479 return ( IsContains( myMeshDS,myShape,anElem,TopAbs_EDGE ) ||
4480 IsContains( myMeshDS,myShape,anElem,TopAbs_FACE ) ||
4481 IsContains( myMeshDS,myShape,anElem,TopAbs_SOLID )||
4482 IsContains( myMeshDS,myShape,anElem,TopAbs_SHELL ));
4484 else if( myType == anElem->GetType() )
4488 case SMDSAbs_Edge : return ( IsContains( myMeshDS,myShape,anElem,TopAbs_EDGE ));
4489 case SMDSAbs_Face : return ( IsContains( myMeshDS,myShape,anElem,TopAbs_FACE ));
4490 case SMDSAbs_Volume: return ( IsContains( myMeshDS,myShape,anElem,TopAbs_SOLID )||
4491 IsContains( myMeshDS,myShape,anElem,TopAbs_SHELL ));
4500 void BelongToGeom::SetType (SMDSAbs_ElementType theType)
4506 SMDSAbs_ElementType BelongToGeom::GetType() const
4511 TopoDS_Shape BelongToGeom::GetShape()
4516 const SMESHDS_Mesh* BelongToGeom::GetMeshDS() const
4521 void BelongToGeom::SetTolerance (double theTolerance)
4523 myTolerance = theTolerance;
4528 double BelongToGeom::GetTolerance()
4535 Description : Predicate for verifying whether entiy lying or partially lying on
4536 specified geometrical support
4539 LyingOnGeom::LyingOnGeom()
4541 myType(SMDSAbs_All),
4542 myIsSubshape(false),
4543 myTolerance(Precision::Confusion())
4546 void LyingOnGeom::SetMesh( const SMDS_Mesh* theMesh )
4548 myMeshDS = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
4552 void LyingOnGeom::SetGeom( const TopoDS_Shape& theShape )
4558 void LyingOnGeom::init()
4560 if (!myMeshDS || myShape.IsNull()) return;
4562 // is sub-shape of main shape?
4563 TopoDS_Shape aMainShape = myMeshDS->ShapeToMesh();
4564 if (aMainShape.IsNull()) {
4565 myIsSubshape = false;
4568 myIsSubshape = myMeshDS->IsGroupOfSubShapes( myShape );
4573 TopTools_IndexedMapOfShape shapes;
4574 TopExp::MapShapes( myShape, shapes );
4575 mySubShapesIDs.Clear();
4576 for ( int i = 1; i <= shapes.Extent(); ++i )
4578 int subID = myMeshDS->ShapeToIndex( shapes( i ));
4580 mySubShapesIDs.Add( subID );
4585 myElementsOnShapePtr.reset(new ElementsOnShape());
4586 myElementsOnShapePtr->SetTolerance(myTolerance);
4587 myElementsOnShapePtr->SetAllNodes(false); // lays on, while true means "belong"
4588 myElementsOnShapePtr->SetMesh(myMeshDS);
4589 myElementsOnShapePtr->SetShape(myShape, myType);
4593 bool LyingOnGeom::IsSatisfy( long theId )
4595 if ( myMeshDS == 0 || myShape.IsNull() )
4600 return myElementsOnShapePtr->IsSatisfy(theId);
4605 const SMDS_MeshElement* elem =
4606 ( myType == SMDSAbs_Node ) ? myMeshDS->FindNode( theId ) : myMeshDS->FindElement( theId );
4608 if ( mySubShapesIDs.Contains( elem->getshapeId() ))
4611 if ( elem->GetType() != SMDSAbs_Node )
4613 SMDS_ElemIteratorPtr nodeItr = elem->nodesIterator();
4614 while ( nodeItr->more() )
4616 const SMDS_MeshElement* aNode = nodeItr->next();
4617 if ( mySubShapesIDs.Contains( aNode->getshapeId() ))
4625 void LyingOnGeom::SetType( SMDSAbs_ElementType theType )
4631 SMDSAbs_ElementType LyingOnGeom::GetType() const
4636 TopoDS_Shape LyingOnGeom::GetShape()
4641 const SMESHDS_Mesh* LyingOnGeom::GetMeshDS() const
4646 void LyingOnGeom::SetTolerance (double theTolerance)
4648 myTolerance = theTolerance;
4653 double LyingOnGeom::GetTolerance()
4658 bool LyingOnGeom::Contains( const SMESHDS_Mesh* theMeshDS,
4659 const TopoDS_Shape& theShape,
4660 const SMDS_MeshElement* theElem,
4661 TopAbs_ShapeEnum theFindShapeEnum,
4662 TopAbs_ShapeEnum theAvoidShapeEnum )
4664 // if (IsContains(theMeshDS, theShape, theElem, theFindShapeEnum, theAvoidShapeEnum))
4667 // TopTools_MapOfShape aSubShapes;
4668 // TopExp_Explorer exp( theShape, theFindShapeEnum, theAvoidShapeEnum );
4669 // for ( ; exp.More(); exp.Next() )
4671 // const TopoDS_Shape& aShape = exp.Current();
4672 // if ( !aSubShapes.Add( aShape )) continue;
4674 // if ( SMESHDS_SubMesh* aSubMesh = theMeshDS->MeshElements( aShape ))
4676 // if ( aSubMesh->Contains( theElem ))
4679 // SMDS_ElemIteratorPtr nodeItr = theElem->nodesIterator();
4680 // while ( nodeItr->more() )
4682 // const SMDS_MeshElement* aNode = nodeItr->next();
4683 // if ( aSubMesh->Contains( aNode ))
4691 TSequenceOfXYZ::TSequenceOfXYZ(): myElem(0)
4694 TSequenceOfXYZ::TSequenceOfXYZ(size_type n) : myArray(n), myElem(0)
4697 TSequenceOfXYZ::TSequenceOfXYZ(size_type n, const gp_XYZ& t) : myArray(n,t), myElem(0)
4700 TSequenceOfXYZ::TSequenceOfXYZ(const TSequenceOfXYZ& theSequenceOfXYZ) : myArray(theSequenceOfXYZ.myArray), myElem(theSequenceOfXYZ.myElem)
4703 template <class InputIterator>
4704 TSequenceOfXYZ::TSequenceOfXYZ(InputIterator theBegin, InputIterator theEnd): myArray(theBegin,theEnd), myElem(0)
4707 TSequenceOfXYZ::~TSequenceOfXYZ()
4710 TSequenceOfXYZ& TSequenceOfXYZ::operator=(const TSequenceOfXYZ& theSequenceOfXYZ)
4712 myArray = theSequenceOfXYZ.myArray;
4713 myElem = theSequenceOfXYZ.myElem;
4717 gp_XYZ& TSequenceOfXYZ::operator()(size_type n)
4719 return myArray[n-1];
4722 const gp_XYZ& TSequenceOfXYZ::operator()(size_type n) const
4724 return myArray[n-1];
4727 void TSequenceOfXYZ::clear()
4732 void TSequenceOfXYZ::reserve(size_type n)
4737 void TSequenceOfXYZ::push_back(const gp_XYZ& v)
4739 myArray.push_back(v);
4742 TSequenceOfXYZ::size_type TSequenceOfXYZ::size() const
4744 return myArray.size();
4747 SMDSAbs_EntityType TSequenceOfXYZ::getElementEntity() const
4749 return myElem ? myElem->GetEntityType() : SMDSEntity_Last;
4752 TMeshModifTracer::TMeshModifTracer():
4753 myMeshModifTime(0), myMesh(0)
4756 void TMeshModifTracer::SetMesh( const SMDS_Mesh* theMesh )
4758 if ( theMesh != myMesh )
4759 myMeshModifTime = 0;
4762 bool TMeshModifTracer::IsMeshModified()
4764 bool modified = false;
4767 modified = ( myMeshModifTime != myMesh->GetMTime() );
4768 myMeshModifTime = myMesh->GetMTime();