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_NodeIteratorPtr anIter = aFace->interlacedNodesIterator();
2436 int i = 0, nbNodes = aFace->NbNodes();
2437 std::vector <const SMDS_MeshNode*> aNodes( nbNodes+1 );
2438 while( anIter->more() )
2439 if ( ! ( aNodes[ i++ ] = anIter->next() ))
2441 aNodes[ nbNodes ] = aNodes[ 0 ];
2443 for ( i = 0; i < nbNodes; i++ )
2444 if ( IsFreeEdge( &aNodes[ i ], theId ) )
2450 SMDSAbs_ElementType FreeEdges::GetType() const
2452 return SMDSAbs_Face;
2455 FreeEdges::Border::Border(long theElemId, long thePntId1, long thePntId2):
2458 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
2459 if(thePntId1 > thePntId2){
2460 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
2464 bool FreeEdges::Border::operator<(const FreeEdges::Border& x) const{
2465 if(myPntId[0] < x.myPntId[0]) return true;
2466 if(myPntId[0] == x.myPntId[0])
2467 if(myPntId[1] < x.myPntId[1]) return true;
2471 inline void UpdateBorders(const FreeEdges::Border& theBorder,
2472 FreeEdges::TBorders& theRegistry,
2473 FreeEdges::TBorders& theContainer)
2475 if(theRegistry.find(theBorder) == theRegistry.end()){
2476 theRegistry.insert(theBorder);
2477 theContainer.insert(theBorder);
2479 theContainer.erase(theBorder);
2483 void FreeEdges::GetBoreders(TBorders& theBorders)
2486 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
2487 for(; anIter->more(); ){
2488 const SMDS_MeshFace* anElem = anIter->next();
2489 long anElemId = anElem->GetID();
2490 SMDS_ElemIteratorPtr aNodesIter;
2491 if ( anElem->IsQuadratic() )
2492 aNodesIter = static_cast<const SMDS_VtkFace*>(anElem)->
2493 interlacedNodesElemIterator();
2495 aNodesIter = anElem->nodesIterator();
2497 const SMDS_MeshElement* aNode;
2498 if(aNodesIter->more()){
2499 aNode = aNodesIter->next();
2500 aNodeId[0] = aNodeId[1] = aNode->GetID();
2502 for(; aNodesIter->more(); ){
2503 aNode = aNodesIter->next();
2504 long anId = aNode->GetID();
2505 Border aBorder(anElemId,aNodeId[1],anId);
2507 UpdateBorders(aBorder,aRegistry,theBorders);
2509 Border aBorder(anElemId,aNodeId[0],aNodeId[1]);
2510 UpdateBorders(aBorder,aRegistry,theBorders);
2514 //================================================================================
2517 Description : Predicate for free nodes
2519 //================================================================================
2521 FreeNodes::FreeNodes()
2526 void FreeNodes::SetMesh( const SMDS_Mesh* theMesh )
2531 bool FreeNodes::IsSatisfy( long theNodeId )
2533 const SMDS_MeshNode* aNode = myMesh->FindNode( theNodeId );
2537 return (aNode->NbInverseElements() < 1);
2540 SMDSAbs_ElementType FreeNodes::GetType() const
2542 return SMDSAbs_Node;
2546 //================================================================================
2549 Description : Predicate for free faces
2551 //================================================================================
2553 FreeFaces::FreeFaces()
2558 void FreeFaces::SetMesh( const SMDS_Mesh* theMesh )
2563 bool FreeFaces::IsSatisfy( long theId )
2565 if (!myMesh) return false;
2566 // check that faces nodes refers to less than two common volumes
2567 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2568 if ( !aFace || aFace->GetType() != SMDSAbs_Face )
2571 int nbNode = aFace->NbNodes();
2573 // collect volumes to check that number of volumes with count equal nbNode not less than 2
2574 typedef map< SMDS_MeshElement*, int > TMapOfVolume; // map of volume counters
2575 typedef map< SMDS_MeshElement*, int >::iterator TItrMapOfVolume; // iterator
2576 TMapOfVolume mapOfVol;
2578 SMDS_ElemIteratorPtr nodeItr = aFace->nodesIterator();
2579 while ( nodeItr->more() ) {
2580 const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>(nodeItr->next());
2581 if ( !aNode ) continue;
2582 SMDS_ElemIteratorPtr volItr = aNode->GetInverseElementIterator(SMDSAbs_Volume);
2583 while ( volItr->more() ) {
2584 SMDS_MeshElement* aVol = (SMDS_MeshElement*)volItr->next();
2585 TItrMapOfVolume itr = mapOfVol.insert(make_pair(aVol, 0)).first;
2590 TItrMapOfVolume volItr = mapOfVol.begin();
2591 TItrMapOfVolume volEnd = mapOfVol.end();
2592 for ( ; volItr != volEnd; ++volItr )
2593 if ( (*volItr).second >= nbNode )
2595 // face is not free if number of volumes constructed on thier nodes more than one
2599 SMDSAbs_ElementType FreeFaces::GetType() const
2601 return SMDSAbs_Face;
2604 //================================================================================
2606 Class : LinearOrQuadratic
2607 Description : Predicate to verify whether a mesh element is linear
2609 //================================================================================
2611 LinearOrQuadratic::LinearOrQuadratic()
2616 void LinearOrQuadratic::SetMesh( const SMDS_Mesh* theMesh )
2621 bool LinearOrQuadratic::IsSatisfy( long theId )
2623 if (!myMesh) return false;
2624 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2625 if ( !anElem || (myType != SMDSAbs_All && anElem->GetType() != myType) )
2627 return (!anElem->IsQuadratic());
2630 void LinearOrQuadratic::SetType( SMDSAbs_ElementType theType )
2635 SMDSAbs_ElementType LinearOrQuadratic::GetType() const
2640 //================================================================================
2643 Description : Functor for check color of group to whic mesh element belongs to
2645 //================================================================================
2647 GroupColor::GroupColor()
2651 bool GroupColor::IsSatisfy( long theId )
2653 return (myIDs.find( theId ) != myIDs.end());
2656 void GroupColor::SetType( SMDSAbs_ElementType theType )
2661 SMDSAbs_ElementType GroupColor::GetType() const
2666 static bool isEqual( const Quantity_Color& theColor1,
2667 const Quantity_Color& theColor2 )
2669 // tolerance to compare colors
2670 const double tol = 5*1e-3;
2671 return ( fabs( theColor1.Red() - theColor2.Red() ) < tol &&
2672 fabs( theColor1.Green() - theColor2.Green() ) < tol &&
2673 fabs( theColor1.Blue() - theColor2.Blue() ) < tol );
2676 void GroupColor::SetMesh( const SMDS_Mesh* theMesh )
2680 const SMESHDS_Mesh* aMesh = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
2684 int nbGrp = aMesh->GetNbGroups();
2688 // iterates on groups and find necessary elements ids
2689 const std::set<SMESHDS_GroupBase*>& aGroups = aMesh->GetGroups();
2690 set<SMESHDS_GroupBase*>::const_iterator GrIt = aGroups.begin();
2691 for (; GrIt != aGroups.end(); GrIt++) {
2692 SMESHDS_GroupBase* aGrp = (*GrIt);
2695 // check type and color of group
2696 if ( !isEqual( myColor, aGrp->GetColor() ) )
2698 if ( myType != SMDSAbs_All && myType != (SMDSAbs_ElementType)aGrp->GetType() )
2701 SMDSAbs_ElementType aGrpElType = (SMDSAbs_ElementType)aGrp->GetType();
2702 if ( myType == aGrpElType || (myType == SMDSAbs_All && aGrpElType != SMDSAbs_Node) ) {
2703 // add elements IDS into control
2704 int aSize = aGrp->Extent();
2705 for (int i = 0; i < aSize; i++)
2706 myIDs.insert( aGrp->GetID(i+1) );
2711 void GroupColor::SetColorStr( const TCollection_AsciiString& theStr )
2713 Kernel_Utils::Localizer loc;
2714 TCollection_AsciiString aStr = theStr;
2715 aStr.RemoveAll( ' ' );
2716 aStr.RemoveAll( '\t' );
2717 for ( int aPos = aStr.Search( ";;" ); aPos != -1; aPos = aStr.Search( ";;" ) )
2718 aStr.Remove( aPos, 2 );
2719 Standard_Real clr[3];
2720 clr[0] = clr[1] = clr[2] = 0.;
2721 for ( int i = 0; i < 3; i++ ) {
2722 TCollection_AsciiString tmpStr = aStr.Token( ";", i+1 );
2723 if ( !tmpStr.IsEmpty() && tmpStr.IsRealValue() )
2724 clr[i] = tmpStr.RealValue();
2726 myColor = Quantity_Color( clr[0], clr[1], clr[2], Quantity_TOC_RGB );
2729 //=======================================================================
2730 // name : GetRangeStr
2731 // Purpose : Get range as a string.
2732 // Example: "1,2,3,50-60,63,67,70-"
2733 //=======================================================================
2735 void GroupColor::GetColorStr( TCollection_AsciiString& theResStr ) const
2738 theResStr += TCollection_AsciiString( myColor.Red() );
2739 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Green() );
2740 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Blue() );
2743 //================================================================================
2745 Class : ElemGeomType
2746 Description : Predicate to check element geometry type
2748 //================================================================================
2750 ElemGeomType::ElemGeomType()
2753 myType = SMDSAbs_All;
2754 myGeomType = SMDSGeom_TRIANGLE;
2757 void ElemGeomType::SetMesh( const SMDS_Mesh* theMesh )
2762 bool ElemGeomType::IsSatisfy( long theId )
2764 if (!myMesh) return false;
2765 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2768 const SMDSAbs_ElementType anElemType = anElem->GetType();
2769 if ( myType != SMDSAbs_All && anElemType != myType )
2771 bool isOk = ( anElem->GetGeomType() == myGeomType );
2775 void ElemGeomType::SetType( SMDSAbs_ElementType theType )
2780 SMDSAbs_ElementType ElemGeomType::GetType() const
2785 void ElemGeomType::SetGeomType( SMDSAbs_GeometryType theType )
2787 myGeomType = theType;
2790 SMDSAbs_GeometryType ElemGeomType::GetGeomType() const
2795 //================================================================================
2797 Class : ElemEntityType
2798 Description : Predicate to check element entity type
2800 //================================================================================
2802 ElemEntityType::ElemEntityType():
2804 myType( SMDSAbs_All ),
2805 myEntityType( SMDSEntity_0D )
2809 void ElemEntityType::SetMesh( const SMDS_Mesh* theMesh )
2814 bool ElemEntityType::IsSatisfy( long theId )
2816 if ( !myMesh ) return false;
2817 if ( myType == SMDSAbs_Node )
2818 return myMesh->FindNode( theId );
2819 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2821 myEntityType == anElem->GetEntityType() );
2824 void ElemEntityType::SetType( SMDSAbs_ElementType theType )
2829 SMDSAbs_ElementType ElemEntityType::GetType() const
2834 void ElemEntityType::SetElemEntityType( SMDSAbs_EntityType theEntityType )
2836 myEntityType = theEntityType;
2839 SMDSAbs_EntityType ElemEntityType::GetElemEntityType() const
2841 return myEntityType;
2844 //================================================================================
2846 * \brief Class ConnectedElements
2848 //================================================================================
2850 ConnectedElements::ConnectedElements():
2851 myNodeID(0), myType( SMDSAbs_All ), myOkIDsReady( false ) {}
2853 SMDSAbs_ElementType ConnectedElements::GetType() const
2856 int ConnectedElements::GetNode() const
2857 { return myXYZ.empty() ? myNodeID : 0; } // myNodeID can be found by myXYZ
2859 std::vector<double> ConnectedElements::GetPoint() const
2862 void ConnectedElements::clearOkIDs()
2863 { myOkIDsReady = false; myOkIDs.clear(); }
2865 void ConnectedElements::SetType( SMDSAbs_ElementType theType )
2867 if ( myType != theType || myMeshModifTracer.IsMeshModified() )
2872 void ConnectedElements::SetMesh( const SMDS_Mesh* theMesh )
2874 myMeshModifTracer.SetMesh( theMesh );
2875 if ( myMeshModifTracer.IsMeshModified() )
2878 if ( !myXYZ.empty() )
2879 SetPoint( myXYZ[0], myXYZ[1], myXYZ[2] ); // find a node near myXYZ it in a new mesh
2883 void ConnectedElements::SetNode( int nodeID )
2888 bool isSameDomain = false;
2889 if ( myOkIDsReady && myMeshModifTracer.GetMesh() && !myMeshModifTracer.IsMeshModified() )
2890 if ( const SMDS_MeshNode* n = myMeshModifTracer.GetMesh()->FindNode( myNodeID ))
2892 SMDS_ElemIteratorPtr eIt = n->GetInverseElementIterator( myType );
2893 while ( !isSameDomain && eIt->more() )
2894 isSameDomain = IsSatisfy( eIt->next()->GetID() );
2896 if ( !isSameDomain )
2900 void ConnectedElements::SetPoint( double x, double y, double z )
2908 bool isSameDomain = false;
2910 // find myNodeID by myXYZ if possible
2911 if ( myMeshModifTracer.GetMesh() )
2913 auto_ptr<SMESH_ElementSearcher> searcher
2914 ( SMESH_MeshAlgos::GetElementSearcher( (SMDS_Mesh&) *myMeshModifTracer.GetMesh() ));
2916 vector< const SMDS_MeshElement* > foundElems;
2917 searcher->FindElementsByPoint( gp_Pnt(x,y,z), SMDSAbs_All, foundElems );
2919 if ( !foundElems.empty() )
2921 myNodeID = foundElems[0]->GetNode(0)->GetID();
2922 if ( myOkIDsReady && !myMeshModifTracer.IsMeshModified() )
2923 isSameDomain = IsSatisfy( foundElems[0]->GetID() );
2926 if ( !isSameDomain )
2930 bool ConnectedElements::IsSatisfy( long theElementId )
2932 // Here we do NOT check if the mesh has changed, we do it in Set...() only!!!
2934 if ( !myOkIDsReady )
2936 if ( !myMeshModifTracer.GetMesh() )
2938 const SMDS_MeshNode* node0 = myMeshModifTracer.GetMesh()->FindNode( myNodeID );
2942 list< const SMDS_MeshNode* > nodeQueue( 1, node0 );
2943 std::set< int > checkedNodeIDs;
2945 // foreach node in nodeQueue:
2946 // foreach element sharing a node:
2947 // add ID of an element of myType to myOkIDs;
2948 // push all element nodes absent from checkedNodeIDs to nodeQueue;
2949 while ( !nodeQueue.empty() )
2951 const SMDS_MeshNode* node = nodeQueue.front();
2952 nodeQueue.pop_front();
2954 // loop on elements sharing the node
2955 SMDS_ElemIteratorPtr eIt = node->GetInverseElementIterator();
2956 while ( eIt->more() )
2958 // keep elements of myType
2959 const SMDS_MeshElement* element = eIt->next();
2960 if ( element->GetType() == myType )
2961 myOkIDs.insert( myOkIDs.end(), element->GetID() );
2963 // enqueue nodes of the element
2964 SMDS_ElemIteratorPtr nIt = element->nodesIterator();
2965 while ( nIt->more() )
2967 const SMDS_MeshNode* n = static_cast< const SMDS_MeshNode* >( nIt->next() );
2968 if ( checkedNodeIDs.insert( n->GetID() ).second )
2969 nodeQueue.push_back( n );
2973 if ( myType == SMDSAbs_Node )
2974 std::swap( myOkIDs, checkedNodeIDs );
2976 size_t totalNbElems = myMeshModifTracer.GetMesh()->GetMeshInfo().NbElements( myType );
2977 if ( myOkIDs.size() == totalNbElems )
2980 myOkIDsReady = true;
2983 return myOkIDs.empty() ? true : myOkIDs.count( theElementId );
2986 //================================================================================
2988 * \brief Class CoplanarFaces
2990 //================================================================================
2992 CoplanarFaces::CoplanarFaces()
2993 : myFaceID(0), myToler(0)
2996 void CoplanarFaces::SetMesh( const SMDS_Mesh* theMesh )
2998 myMeshModifTracer.SetMesh( theMesh );
2999 if ( myMeshModifTracer.IsMeshModified() )
3001 // Build a set of coplanar face ids
3003 myCoplanarIDs.clear();
3005 if ( !myMeshModifTracer.GetMesh() || !myFaceID || !myToler )
3008 const SMDS_MeshElement* face = myMeshModifTracer.GetMesh()->FindElement( myFaceID );
3009 if ( !face || face->GetType() != SMDSAbs_Face )
3013 gp_Vec myNorm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
3017 const double radianTol = myToler * M_PI / 180.;
3018 std::set< SMESH_TLink > checkedLinks;
3020 std::list< pair< const SMDS_MeshElement*, gp_Vec > > faceQueue;
3021 faceQueue.push_back( make_pair( face, myNorm ));
3022 while ( !faceQueue.empty() )
3024 face = faceQueue.front().first;
3025 myNorm = faceQueue.front().second;
3026 faceQueue.pop_front();
3028 for ( int i = 0, nbN = face->NbCornerNodes(); i < nbN; ++i )
3030 const SMDS_MeshNode* n1 = face->GetNode( i );
3031 const SMDS_MeshNode* n2 = face->GetNode(( i+1 )%nbN);
3032 if ( !checkedLinks.insert( SMESH_TLink( n1, n2 )).second )
3034 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator(SMDSAbs_Face);
3035 while ( fIt->more() )
3037 const SMDS_MeshElement* f = fIt->next();
3038 if ( f->GetNodeIndex( n2 ) > -1 )
3040 gp_Vec norm = getNormale( static_cast<const SMDS_MeshFace*>(f), &normOK );
3041 if (!normOK || myNorm.Angle( norm ) <= radianTol)
3043 myCoplanarIDs.insert( f->GetID() );
3044 faceQueue.push_back( make_pair( f, norm ));
3052 bool CoplanarFaces::IsSatisfy( long theElementId )
3054 return myCoplanarIDs.count( theElementId );
3059 *Description : Predicate for Range of Ids.
3060 * Range may be specified with two ways.
3061 * 1. Using AddToRange method
3062 * 2. With SetRangeStr method. Parameter of this method is a string
3063 * like as "1,2,3,50-60,63,67,70-"
3066 //=======================================================================
3067 // name : RangeOfIds
3068 // Purpose : Constructor
3069 //=======================================================================
3070 RangeOfIds::RangeOfIds()
3073 myType = SMDSAbs_All;
3076 //=======================================================================
3078 // Purpose : Set mesh
3079 //=======================================================================
3080 void RangeOfIds::SetMesh( const SMDS_Mesh* theMesh )
3085 //=======================================================================
3086 // name : AddToRange
3087 // Purpose : Add ID to the range
3088 //=======================================================================
3089 bool RangeOfIds::AddToRange( long theEntityId )
3091 myIds.Add( theEntityId );
3095 //=======================================================================
3096 // name : GetRangeStr
3097 // Purpose : Get range as a string.
3098 // Example: "1,2,3,50-60,63,67,70-"
3099 //=======================================================================
3100 void RangeOfIds::GetRangeStr( TCollection_AsciiString& theResStr )
3104 TColStd_SequenceOfInteger anIntSeq;
3105 TColStd_SequenceOfAsciiString aStrSeq;
3107 TColStd_MapIteratorOfMapOfInteger anIter( myIds );
3108 for ( ; anIter.More(); anIter.Next() )
3110 int anId = anIter.Key();
3111 TCollection_AsciiString aStr( anId );
3112 anIntSeq.Append( anId );
3113 aStrSeq.Append( aStr );
3116 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
3118 int aMinId = myMin( i );
3119 int aMaxId = myMax( i );
3121 TCollection_AsciiString aStr;
3122 if ( aMinId != IntegerFirst() )
3127 if ( aMaxId != IntegerLast() )
3130 // find position of the string in result sequence and insert string in it
3131 if ( anIntSeq.Length() == 0 )
3133 anIntSeq.Append( aMinId );
3134 aStrSeq.Append( aStr );
3138 if ( aMinId < anIntSeq.First() )
3140 anIntSeq.Prepend( aMinId );
3141 aStrSeq.Prepend( aStr );
3143 else if ( aMinId > anIntSeq.Last() )
3145 anIntSeq.Append( aMinId );
3146 aStrSeq.Append( aStr );
3149 for ( int j = 1, k = anIntSeq.Length(); j <= k; j++ )
3150 if ( aMinId < anIntSeq( j ) )
3152 anIntSeq.InsertBefore( j, aMinId );
3153 aStrSeq.InsertBefore( j, aStr );
3159 if ( aStrSeq.Length() == 0 )
3162 theResStr = aStrSeq( 1 );
3163 for ( int j = 2, k = aStrSeq.Length(); j <= k; j++ )
3166 theResStr += aStrSeq( j );
3170 //=======================================================================
3171 // name : SetRangeStr
3172 // Purpose : Define range with string
3173 // Example of entry string: "1,2,3,50-60,63,67,70-"
3174 //=======================================================================
3175 bool RangeOfIds::SetRangeStr( const TCollection_AsciiString& theStr )
3181 TCollection_AsciiString aStr = theStr;
3182 //aStr.RemoveAll( ' ' );
3183 //aStr.RemoveAll( '\t' );
3184 for ( int i = 1; i <= aStr.Length(); ++i )
3185 if ( isspace( aStr.Value( i )))
3186 aStr.SetValue( i, ',');
3188 for ( int aPos = aStr.Search( ",," ); aPos != -1; aPos = aStr.Search( ",," ) )
3189 aStr.Remove( aPos, 1 );
3191 TCollection_AsciiString tmpStr = aStr.Token( ",", 1 );
3193 while ( tmpStr != "" )
3195 tmpStr = aStr.Token( ",", i++ );
3196 int aPos = tmpStr.Search( '-' );
3200 if ( tmpStr.IsIntegerValue() )
3201 myIds.Add( tmpStr.IntegerValue() );
3207 TCollection_AsciiString aMaxStr = tmpStr.Split( aPos );
3208 TCollection_AsciiString aMinStr = tmpStr;
3210 while ( aMinStr.Search( "-" ) != -1 ) aMinStr.RemoveAll( '-' );
3211 while ( aMaxStr.Search( "-" ) != -1 ) aMaxStr.RemoveAll( '-' );
3213 if ( (!aMinStr.IsEmpty() && !aMinStr.IsIntegerValue()) ||
3214 (!aMaxStr.IsEmpty() && !aMaxStr.IsIntegerValue()) )
3217 myMin.Append( aMinStr.IsEmpty() ? IntegerFirst() : aMinStr.IntegerValue() );
3218 myMax.Append( aMaxStr.IsEmpty() ? IntegerLast() : aMaxStr.IntegerValue() );
3225 //=======================================================================
3227 // Purpose : Get type of supported entities
3228 //=======================================================================
3229 SMDSAbs_ElementType RangeOfIds::GetType() const
3234 //=======================================================================
3236 // Purpose : Set type of supported entities
3237 //=======================================================================
3238 void RangeOfIds::SetType( SMDSAbs_ElementType theType )
3243 //=======================================================================
3245 // Purpose : Verify whether entity satisfies to this rpedicate
3246 //=======================================================================
3247 bool RangeOfIds::IsSatisfy( long theId )
3252 if ( myType == SMDSAbs_Node )
3254 if ( myMesh->FindNode( theId ) == 0 )
3259 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
3260 if ( anElem == 0 || (myType != anElem->GetType() && myType != SMDSAbs_All ))
3264 if ( myIds.Contains( theId ) )
3267 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
3268 if ( theId >= myMin( i ) && theId <= myMax( i ) )
3276 Description : Base class for comparators
3278 Comparator::Comparator():
3282 Comparator::~Comparator()
3285 void Comparator::SetMesh( const SMDS_Mesh* theMesh )
3288 myFunctor->SetMesh( theMesh );
3291 void Comparator::SetMargin( double theValue )
3293 myMargin = theValue;
3296 void Comparator::SetNumFunctor( NumericalFunctorPtr theFunct )
3298 myFunctor = theFunct;
3301 SMDSAbs_ElementType Comparator::GetType() const
3303 return myFunctor ? myFunctor->GetType() : SMDSAbs_All;
3306 double Comparator::GetMargin()
3314 Description : Comparator "<"
3316 bool LessThan::IsSatisfy( long theId )
3318 return myFunctor && myFunctor->GetValue( theId ) < myMargin;
3324 Description : Comparator ">"
3326 bool MoreThan::IsSatisfy( long theId )
3328 return myFunctor && myFunctor->GetValue( theId ) > myMargin;
3334 Description : Comparator "="
3337 myToler(Precision::Confusion())
3340 bool EqualTo::IsSatisfy( long theId )
3342 return myFunctor && fabs( myFunctor->GetValue( theId ) - myMargin ) < myToler;
3345 void EqualTo::SetTolerance( double theToler )
3350 double EqualTo::GetTolerance()
3357 Description : Logical NOT predicate
3359 LogicalNOT::LogicalNOT()
3362 LogicalNOT::~LogicalNOT()
3365 bool LogicalNOT::IsSatisfy( long theId )
3367 return myPredicate && !myPredicate->IsSatisfy( theId );
3370 void LogicalNOT::SetMesh( const SMDS_Mesh* theMesh )
3373 myPredicate->SetMesh( theMesh );
3376 void LogicalNOT::SetPredicate( PredicatePtr thePred )
3378 myPredicate = thePred;
3381 SMDSAbs_ElementType LogicalNOT::GetType() const
3383 return myPredicate ? myPredicate->GetType() : SMDSAbs_All;
3388 Class : LogicalBinary
3389 Description : Base class for binary logical predicate
3391 LogicalBinary::LogicalBinary()
3394 LogicalBinary::~LogicalBinary()
3397 void LogicalBinary::SetMesh( const SMDS_Mesh* theMesh )
3400 myPredicate1->SetMesh( theMesh );
3403 myPredicate2->SetMesh( theMesh );
3406 void LogicalBinary::SetPredicate1( PredicatePtr thePredicate )
3408 myPredicate1 = thePredicate;
3411 void LogicalBinary::SetPredicate2( PredicatePtr thePredicate )
3413 myPredicate2 = thePredicate;
3416 SMDSAbs_ElementType LogicalBinary::GetType() const
3418 if ( !myPredicate1 || !myPredicate2 )
3421 SMDSAbs_ElementType aType1 = myPredicate1->GetType();
3422 SMDSAbs_ElementType aType2 = myPredicate2->GetType();
3424 return aType1 == aType2 ? aType1 : SMDSAbs_All;
3430 Description : Logical AND
3432 bool LogicalAND::IsSatisfy( long theId )
3437 myPredicate1->IsSatisfy( theId ) &&
3438 myPredicate2->IsSatisfy( theId );
3444 Description : Logical OR
3446 bool LogicalOR::IsSatisfy( long theId )
3451 (myPredicate1->IsSatisfy( theId ) ||
3452 myPredicate2->IsSatisfy( theId ));
3461 // #include <tbb/parallel_for.h>
3462 // #include <tbb/enumerable_thread_specific.h>
3464 // namespace Parallel
3466 // typedef tbb::enumerable_thread_specific< TIdSequence > TIdSeq;
3470 // const SMDS_Mesh* myMesh;
3471 // PredicatePtr myPredicate;
3472 // TIdSeq & myOKIds;
3473 // Predicate( const SMDS_Mesh* m, PredicatePtr p, TIdSeq & ids ):
3474 // myMesh(m), myPredicate(p->Duplicate()), myOKIds(ids) {}
3475 // void operator() ( const tbb::blocked_range<size_t>& r ) const
3477 // for ( size_t i = r.begin(); i != r.end(); ++i )
3478 // if ( myPredicate->IsSatisfy( i ))
3479 // myOKIds.local().push_back();
3491 void Filter::SetPredicate( PredicatePtr thePredicate )
3493 myPredicate = thePredicate;
3496 void Filter::GetElementsId( const SMDS_Mesh* theMesh,
3497 PredicatePtr thePredicate,
3498 TIdSequence& theSequence )
3500 theSequence.clear();
3502 if ( !theMesh || !thePredicate )
3505 thePredicate->SetMesh( theMesh );
3507 SMDS_ElemIteratorPtr elemIt = theMesh->elementsIterator( thePredicate->GetType() );
3509 while ( elemIt->more() ) {
3510 const SMDS_MeshElement* anElem = elemIt->next();
3511 long anId = anElem->GetID();
3512 if ( thePredicate->IsSatisfy( anId ) )
3513 theSequence.push_back( anId );
3518 void Filter::GetElementsId( const SMDS_Mesh* theMesh,
3519 Filter::TIdSequence& theSequence )
3521 GetElementsId(theMesh,myPredicate,theSequence);
3528 typedef std::set<SMDS_MeshFace*> TMapOfFacePtr;
3534 ManifoldPart::Link::Link( SMDS_MeshNode* theNode1,
3535 SMDS_MeshNode* theNode2 )
3541 ManifoldPart::Link::~Link()
3547 bool ManifoldPart::Link::IsEqual( const ManifoldPart::Link& theLink ) const
3549 if ( myNode1 == theLink.myNode1 &&
3550 myNode2 == theLink.myNode2 )
3552 else if ( myNode1 == theLink.myNode2 &&
3553 myNode2 == theLink.myNode1 )
3559 bool ManifoldPart::Link::operator<( const ManifoldPart::Link& x ) const
3561 if(myNode1 < x.myNode1) return true;
3562 if(myNode1 == x.myNode1)
3563 if(myNode2 < x.myNode2) return true;
3567 bool ManifoldPart::IsEqual( const ManifoldPart::Link& theLink1,
3568 const ManifoldPart::Link& theLink2 )
3570 return theLink1.IsEqual( theLink2 );
3573 ManifoldPart::ManifoldPart()
3576 myAngToler = Precision::Angular();
3577 myIsOnlyManifold = true;
3580 ManifoldPart::~ManifoldPart()
3585 void ManifoldPart::SetMesh( const SMDS_Mesh* theMesh )
3591 SMDSAbs_ElementType ManifoldPart::GetType() const
3592 { return SMDSAbs_Face; }
3594 bool ManifoldPart::IsSatisfy( long theElementId )
3596 return myMapIds.Contains( theElementId );
3599 void ManifoldPart::SetAngleTolerance( const double theAngToler )
3600 { myAngToler = theAngToler; }
3602 double ManifoldPart::GetAngleTolerance() const
3603 { return myAngToler; }
3605 void ManifoldPart::SetIsOnlyManifold( const bool theIsOnly )
3606 { myIsOnlyManifold = theIsOnly; }
3608 void ManifoldPart::SetStartElem( const long theStartId )
3609 { myStartElemId = theStartId; }
3611 bool ManifoldPart::process()
3614 myMapBadGeomIds.Clear();
3616 myAllFacePtr.clear();
3617 myAllFacePtrIntDMap.clear();
3621 // collect all faces into own map
3622 SMDS_FaceIteratorPtr anFaceItr = myMesh->facesIterator();
3623 for (; anFaceItr->more(); )
3625 SMDS_MeshFace* aFacePtr = (SMDS_MeshFace*)anFaceItr->next();
3626 myAllFacePtr.push_back( aFacePtr );
3627 myAllFacePtrIntDMap[aFacePtr] = myAllFacePtr.size()-1;
3630 SMDS_MeshFace* aStartFace = (SMDS_MeshFace*)myMesh->FindElement( myStartElemId );
3634 // the map of non manifold links and bad geometry
3635 TMapOfLink aMapOfNonManifold;
3636 TColStd_MapOfInteger aMapOfTreated;
3638 // begin cycle on faces from start index and run on vector till the end
3639 // and from begin to start index to cover whole vector
3640 const int aStartIndx = myAllFacePtrIntDMap[aStartFace];
3641 bool isStartTreat = false;
3642 for ( int fi = aStartIndx; !isStartTreat || fi != aStartIndx ; fi++ )
3644 if ( fi == aStartIndx )
3645 isStartTreat = true;
3646 // as result next time when fi will be equal to aStartIndx
3648 SMDS_MeshFace* aFacePtr = myAllFacePtr[ fi ];
3649 if ( aMapOfTreated.Contains( aFacePtr->GetID() ) )
3652 aMapOfTreated.Add( aFacePtr->GetID() );
3653 TColStd_MapOfInteger aResFaces;
3654 if ( !findConnected( myAllFacePtrIntDMap, aFacePtr,
3655 aMapOfNonManifold, aResFaces ) )
3657 TColStd_MapIteratorOfMapOfInteger anItr( aResFaces );
3658 for ( ; anItr.More(); anItr.Next() )
3660 int aFaceId = anItr.Key();
3661 aMapOfTreated.Add( aFaceId );
3662 myMapIds.Add( aFaceId );
3665 if ( fi == ( myAllFacePtr.size() - 1 ) )
3667 } // end run on vector of faces
3668 return !myMapIds.IsEmpty();
3671 static void getLinks( const SMDS_MeshFace* theFace,
3672 ManifoldPart::TVectorOfLink& theLinks )
3674 int aNbNode = theFace->NbNodes();
3675 SMDS_ElemIteratorPtr aNodeItr = theFace->nodesIterator();
3677 SMDS_MeshNode* aNode = 0;
3678 for ( ; aNodeItr->more() && i <= aNbNode; )
3681 SMDS_MeshNode* aN1 = (SMDS_MeshNode*)aNodeItr->next();
3685 SMDS_MeshNode* aN2 = ( i >= aNbNode ) ? aNode : (SMDS_MeshNode*)aNodeItr->next();
3687 ManifoldPart::Link aLink( aN1, aN2 );
3688 theLinks.push_back( aLink );
3692 bool ManifoldPart::findConnected
3693 ( const ManifoldPart::TDataMapFacePtrInt& theAllFacePtrInt,
3694 SMDS_MeshFace* theStartFace,
3695 ManifoldPart::TMapOfLink& theNonManifold,
3696 TColStd_MapOfInteger& theResFaces )
3698 theResFaces.Clear();
3699 if ( !theAllFacePtrInt.size() )
3702 if ( getNormale( theStartFace ).SquareModulus() <= gp::Resolution() )
3704 myMapBadGeomIds.Add( theStartFace->GetID() );
3708 ManifoldPart::TMapOfLink aMapOfBoundary, aMapToSkip;
3709 ManifoldPart::TVectorOfLink aSeqOfBoundary;
3710 theResFaces.Add( theStartFace->GetID() );
3711 ManifoldPart::TDataMapOfLinkFacePtr aDMapLinkFace;
3713 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3714 aDMapLinkFace, theNonManifold, theStartFace );
3716 bool isDone = false;
3717 while ( !isDone && aMapOfBoundary.size() != 0 )
3719 bool isToReset = false;
3720 ManifoldPart::TVectorOfLink::iterator pLink = aSeqOfBoundary.begin();
3721 for ( ; !isToReset && pLink != aSeqOfBoundary.end(); ++pLink )
3723 ManifoldPart::Link aLink = *pLink;
3724 if ( aMapToSkip.find( aLink ) != aMapToSkip.end() )
3726 // each link could be treated only once
3727 aMapToSkip.insert( aLink );
3729 ManifoldPart::TVectorOfFacePtr aFaces;
3731 if ( myIsOnlyManifold &&
3732 (theNonManifold.find( aLink ) != theNonManifold.end()) )
3736 getFacesByLink( aLink, aFaces );
3737 // filter the element to keep only indicated elements
3738 ManifoldPart::TVectorOfFacePtr aFiltered;
3739 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3740 for ( ; pFace != aFaces.end(); ++pFace )
3742 SMDS_MeshFace* aFace = *pFace;
3743 if ( myAllFacePtrIntDMap.find( aFace ) != myAllFacePtrIntDMap.end() )
3744 aFiltered.push_back( aFace );
3747 if ( aFaces.size() < 2 ) // no neihgbour faces
3749 else if ( myIsOnlyManifold && aFaces.size() > 2 ) // non manifold case
3751 theNonManifold.insert( aLink );
3756 // compare normal with normals of neighbor element
3757 SMDS_MeshFace* aPrevFace = aDMapLinkFace[ aLink ];
3758 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3759 for ( ; pFace != aFaces.end(); ++pFace )
3761 SMDS_MeshFace* aNextFace = *pFace;
3762 if ( aPrevFace == aNextFace )
3764 int anNextFaceID = aNextFace->GetID();
3765 if ( myIsOnlyManifold && theResFaces.Contains( anNextFaceID ) )
3766 // should not be with non manifold restriction. probably bad topology
3768 // check if face was treated and skipped
3769 if ( myMapBadGeomIds.Contains( anNextFaceID ) ||
3770 !isInPlane( aPrevFace, aNextFace ) )
3772 // add new element to connected and extend the boundaries.
3773 theResFaces.Add( anNextFaceID );
3774 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3775 aDMapLinkFace, theNonManifold, aNextFace );
3779 isDone = !isToReset;
3782 return !theResFaces.IsEmpty();
3785 bool ManifoldPart::isInPlane( const SMDS_MeshFace* theFace1,
3786 const SMDS_MeshFace* theFace2 )
3788 gp_Dir aNorm1 = gp_Dir( getNormale( theFace1 ) );
3789 gp_XYZ aNorm2XYZ = getNormale( theFace2 );
3790 if ( aNorm2XYZ.SquareModulus() <= gp::Resolution() )
3792 myMapBadGeomIds.Add( theFace2->GetID() );
3795 if ( aNorm1.IsParallel( gp_Dir( aNorm2XYZ ), myAngToler ) )
3801 void ManifoldPart::expandBoundary
3802 ( ManifoldPart::TMapOfLink& theMapOfBoundary,
3803 ManifoldPart::TVectorOfLink& theSeqOfBoundary,
3804 ManifoldPart::TDataMapOfLinkFacePtr& theDMapLinkFacePtr,
3805 ManifoldPart::TMapOfLink& theNonManifold,
3806 SMDS_MeshFace* theNextFace ) const
3808 ManifoldPart::TVectorOfLink aLinks;
3809 getLinks( theNextFace, aLinks );
3810 int aNbLink = (int)aLinks.size();
3811 for ( int i = 0; i < aNbLink; i++ )
3813 ManifoldPart::Link aLink = aLinks[ i ];
3814 if ( myIsOnlyManifold && (theNonManifold.find( aLink ) != theNonManifold.end()) )
3816 if ( theMapOfBoundary.find( aLink ) != theMapOfBoundary.end() )
3818 if ( myIsOnlyManifold )
3820 // remove from boundary
3821 theMapOfBoundary.erase( aLink );
3822 ManifoldPart::TVectorOfLink::iterator pLink = theSeqOfBoundary.begin();
3823 for ( ; pLink != theSeqOfBoundary.end(); ++pLink )
3825 ManifoldPart::Link aBoundLink = *pLink;
3826 if ( aBoundLink.IsEqual( aLink ) )
3828 theSeqOfBoundary.erase( pLink );
3836 theMapOfBoundary.insert( aLink );
3837 theSeqOfBoundary.push_back( aLink );
3838 theDMapLinkFacePtr[ aLink ] = theNextFace;
3843 void ManifoldPart::getFacesByLink( const ManifoldPart::Link& theLink,
3844 ManifoldPart::TVectorOfFacePtr& theFaces ) const
3846 std::set<SMDS_MeshCell *> aSetOfFaces;
3847 // take all faces that shared first node
3848 SMDS_ElemIteratorPtr anItr = theLink.myNode1->facesIterator();
3849 for ( ; anItr->more(); )
3851 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3854 aSetOfFaces.insert( aFace );
3856 // take all faces that shared second node
3857 anItr = theLink.myNode2->facesIterator();
3858 // find the common part of two sets
3859 for ( ; anItr->more(); )
3861 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3862 if ( aSetOfFaces.count( aFace ) )
3863 theFaces.push_back( aFace );
3868 Class : BelongToMeshGroup
3869 Description : Verify whether a mesh element is included into a mesh group
3871 BelongToMeshGroup::BelongToMeshGroup(): myGroup( 0 )
3875 void BelongToMeshGroup::SetGroup( SMESHDS_GroupBase* g )
3880 void BelongToMeshGroup::SetStoreName( const std::string& sn )
3885 void BelongToMeshGroup::SetMesh( const SMDS_Mesh* theMesh )
3887 if ( myGroup && myGroup->GetMesh() != theMesh )
3891 if ( !myGroup && !myStoreName.empty() )
3893 if ( const SMESHDS_Mesh* aMesh = dynamic_cast<const SMESHDS_Mesh*>(theMesh))
3895 const std::set<SMESHDS_GroupBase*>& grps = aMesh->GetGroups();
3896 std::set<SMESHDS_GroupBase*>::const_iterator g = grps.begin();
3897 for ( ; g != grps.end() && !myGroup; ++g )
3898 if ( *g && myStoreName == (*g)->GetStoreName() )
3904 myGroup->IsEmpty(); // make GroupOnFilter update its predicate
3908 bool BelongToMeshGroup::IsSatisfy( long theElementId )
3910 return myGroup ? myGroup->Contains( theElementId ) : false;
3913 SMDSAbs_ElementType BelongToMeshGroup::GetType() const
3915 return myGroup ? myGroup->GetType() : SMDSAbs_All;
3922 ElementsOnSurface::ElementsOnSurface()
3925 myType = SMDSAbs_All;
3927 myToler = Precision::Confusion();
3928 myUseBoundaries = false;
3931 ElementsOnSurface::~ElementsOnSurface()
3935 void ElementsOnSurface::SetMesh( const SMDS_Mesh* theMesh )
3937 myMeshModifTracer.SetMesh( theMesh );
3938 if ( myMeshModifTracer.IsMeshModified())
3942 bool ElementsOnSurface::IsSatisfy( long theElementId )
3944 return myIds.Contains( theElementId );
3947 SMDSAbs_ElementType ElementsOnSurface::GetType() const
3950 void ElementsOnSurface::SetTolerance( const double theToler )
3952 if ( myToler != theToler )
3957 double ElementsOnSurface::GetTolerance() const
3960 void ElementsOnSurface::SetUseBoundaries( bool theUse )
3962 if ( myUseBoundaries != theUse ) {
3963 myUseBoundaries = theUse;
3964 SetSurface( mySurf, myType );
3968 void ElementsOnSurface::SetSurface( const TopoDS_Shape& theShape,
3969 const SMDSAbs_ElementType theType )
3974 if ( theShape.IsNull() || theShape.ShapeType() != TopAbs_FACE )
3976 mySurf = TopoDS::Face( theShape );
3977 BRepAdaptor_Surface SA( mySurf, myUseBoundaries );
3979 u1 = SA.FirstUParameter(),
3980 u2 = SA.LastUParameter(),
3981 v1 = SA.FirstVParameter(),
3982 v2 = SA.LastVParameter();
3983 Handle(Geom_Surface) surf = BRep_Tool::Surface( mySurf );
3984 myProjector.Init( surf, u1,u2, v1,v2 );
3988 void ElementsOnSurface::process()
3991 if ( mySurf.IsNull() )
3994 if ( !myMeshModifTracer.GetMesh() )
3997 myIds.ReSize( myMeshModifTracer.GetMesh()->GetMeshInfo().NbElements( myType ));
3999 SMDS_ElemIteratorPtr anIter = myMeshModifTracer.GetMesh()->elementsIterator( myType );
4000 for(; anIter->more(); )
4001 process( anIter->next() );
4004 void ElementsOnSurface::process( const SMDS_MeshElement* theElemPtr )
4006 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
4007 bool isSatisfy = true;
4008 for ( ; aNodeItr->more(); )
4010 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
4011 if ( !isOnSurface( aNode ) )
4018 myIds.Add( theElemPtr->GetID() );
4021 bool ElementsOnSurface::isOnSurface( const SMDS_MeshNode* theNode )
4023 if ( mySurf.IsNull() )
4026 gp_Pnt aPnt( theNode->X(), theNode->Y(), theNode->Z() );
4027 // double aToler2 = myToler * myToler;
4028 // if ( mySurf->IsKind(STANDARD_TYPE(Geom_Plane)))
4030 // gp_Pln aPln = Handle(Geom_Plane)::DownCast(mySurf)->Pln();
4031 // if ( aPln.SquareDistance( aPnt ) > aToler2 )
4034 // else if ( mySurf->IsKind(STANDARD_TYPE(Geom_CylindricalSurface)))
4036 // gp_Cylinder aCyl = Handle(Geom_CylindricalSurface)::DownCast(mySurf)->Cylinder();
4037 // double aRad = aCyl.Radius();
4038 // gp_Ax3 anAxis = aCyl.Position();
4039 // gp_XYZ aLoc = aCyl.Location().XYZ();
4040 // double aXDist = anAxis.XDirection().XYZ() * ( aPnt.XYZ() - aLoc );
4041 // double aYDist = anAxis.YDirection().XYZ() * ( aPnt.XYZ() - aLoc );
4042 // if ( fabs(aXDist*aXDist + aYDist*aYDist - aRad*aRad) > aToler2 )
4047 myProjector.Perform( aPnt );
4048 bool isOn = ( myProjector.IsDone() && myProjector.LowerDistance() <= myToler );
4058 ElementsOnShape::ElementsOnShape()
4060 myType(SMDSAbs_All),
4061 myToler(Precision::Confusion()),
4062 myAllNodesFlag(false)
4066 ElementsOnShape::~ElementsOnShape()
4071 SMDSAbs_ElementType ElementsOnShape::GetType() const
4076 void ElementsOnShape::SetTolerance (const double theToler)
4078 if (myToler != theToler) {
4080 SetShape(myShape, myType);
4084 double ElementsOnShape::GetTolerance() const
4089 void ElementsOnShape::SetAllNodes (bool theAllNodes)
4091 myAllNodesFlag = theAllNodes;
4094 void ElementsOnShape::SetMesh (const SMDS_Mesh* theMesh)
4096 myMeshModifTracer.SetMesh( theMesh );
4097 if ( myMeshModifTracer.IsMeshModified())
4099 size_t nbNodes = theMesh ? theMesh->NbNodes() : 0;
4100 if ( myNodeIsChecked.size() == nbNodes )
4102 std::fill( myNodeIsChecked.begin(), myNodeIsChecked.end(), false );
4106 SMESHUtils::FreeVector( myNodeIsChecked );
4107 SMESHUtils::FreeVector( myNodeIsOut );
4108 myNodeIsChecked.resize( nbNodes, false );
4109 myNodeIsOut.resize( nbNodes );
4114 bool ElementsOnShape::getNodeIsOut( const SMDS_MeshNode* n, bool& isOut )
4116 if ( n->GetID() >= (int) myNodeIsChecked.size() ||
4117 !myNodeIsChecked[ n->GetID() ])
4120 isOut = myNodeIsOut[ n->GetID() ];
4124 void ElementsOnShape::setNodeIsOut( const SMDS_MeshNode* n, bool isOut )
4126 if ( n->GetID() < (int) myNodeIsChecked.size() )
4128 myNodeIsChecked[ n->GetID() ] = true;
4129 myNodeIsOut [ n->GetID() ] = isOut;
4133 void ElementsOnShape::SetShape (const TopoDS_Shape& theShape,
4134 const SMDSAbs_ElementType theType)
4138 if ( myShape.IsNull() ) return;
4140 TopTools_IndexedMapOfShape shapesMap;
4141 TopAbs_ShapeEnum shapeTypes[4] = { TopAbs_SOLID, TopAbs_FACE, TopAbs_EDGE, TopAbs_VERTEX };
4142 TopExp_Explorer sub;
4143 for ( int i = 0; i < 4; ++i )
4145 if ( shapesMap.IsEmpty() )
4146 for ( sub.Init( myShape, shapeTypes[i] ); sub.More(); sub.Next() )
4147 shapesMap.Add( sub.Current() );
4149 for ( sub.Init( myShape, shapeTypes[i], shapeTypes[i-1] ); sub.More(); sub.Next() )
4150 shapesMap.Add( sub.Current() );
4154 myClassifiers.resize( shapesMap.Extent() );
4155 for ( int i = 0; i < shapesMap.Extent(); ++i )
4156 myClassifiers[ i ] = new TClassifier( shapesMap( i+1 ), myToler );
4158 if ( theType == SMDSAbs_Node )
4160 SMESHUtils::FreeVector( myNodeIsChecked );
4161 SMESHUtils::FreeVector( myNodeIsOut );
4165 std::fill( myNodeIsChecked.begin(), myNodeIsChecked.end(), false );
4169 void ElementsOnShape::clearClassifiers()
4171 for ( size_t i = 0; i < myClassifiers.size(); ++i )
4172 delete myClassifiers[ i ];
4173 myClassifiers.clear();
4176 bool ElementsOnShape::IsSatisfy (long elemId)
4178 const SMDS_Mesh* mesh = myMeshModifTracer.GetMesh();
4179 const SMDS_MeshElement* elem =
4180 ( myType == SMDSAbs_Node ? mesh->FindNode( elemId ) : mesh->FindElement( elemId ));
4181 if ( !elem || myClassifiers.empty() )
4184 bool isSatisfy = myAllNodesFlag, isNodeOut;
4186 gp_XYZ centerXYZ (0, 0, 0);
4188 SMDS_ElemIteratorPtr aNodeItr = elem->nodesIterator();
4189 while (aNodeItr->more() && (isSatisfy == myAllNodesFlag))
4191 SMESH_TNodeXYZ aPnt( aNodeItr->next() );
4195 if ( !getNodeIsOut( aPnt._node, isNodeOut ))
4197 for ( size_t i = 0; i < myClassifiers.size() && isNodeOut; ++i )
4198 isNodeOut = myClassifiers[i]->IsOut( aPnt );
4200 setNodeIsOut( aPnt._node, isNodeOut );
4202 isSatisfy = !isNodeOut;
4205 // Check the center point for volumes MantisBug 0020168
4208 myClassifiers[0]->ShapeType() == TopAbs_SOLID)
4210 centerXYZ /= elem->NbNodes();
4212 for ( size_t i = 0; i < myClassifiers.size() && !isSatisfy; ++i )
4213 isSatisfy = ! myClassifiers[i]->IsOut( centerXYZ );
4219 TopAbs_ShapeEnum ElementsOnShape::TClassifier::ShapeType() const
4221 return myShape.ShapeType();
4224 bool ElementsOnShape::TClassifier::IsOut(const gp_Pnt& p)
4226 return (this->*myIsOutFun)( p );
4229 void ElementsOnShape::TClassifier::Init (const TopoDS_Shape& theShape, double theTol)
4233 switch ( myShape.ShapeType() )
4235 case TopAbs_SOLID: {
4236 if ( isBox( theShape ))
4238 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfBox;
4242 mySolidClfr.Load(theShape);
4243 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfSolid;
4248 Standard_Real u1,u2,v1,v2;
4249 Handle(Geom_Surface) surf = BRep_Tool::Surface( TopoDS::Face( theShape ));
4250 surf->Bounds( u1,u2,v1,v2 );
4251 myProjFace.Init(surf, u1,u2, v1,v2, myTol );
4252 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfFace;
4256 Standard_Real u1, u2;
4257 Handle(Geom_Curve) curve = BRep_Tool::Curve( TopoDS::Edge(theShape), u1, u2);
4258 myProjEdge.Init(curve, u1, u2);
4259 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfEdge;
4262 case TopAbs_VERTEX:{
4263 myVertexXYZ = BRep_Tool::Pnt( TopoDS::Vertex( theShape ) );
4264 myIsOutFun = & ElementsOnShape::TClassifier::isOutOfVertex;
4268 throw SALOME_Exception("Programmer error in usage of ElementsOnShape::TClassifier");
4272 bool ElementsOnShape::TClassifier::isOutOfSolid (const gp_Pnt& p)
4274 mySolidClfr.Perform( p, myTol );
4275 return ( mySolidClfr.State() != TopAbs_IN && mySolidClfr.State() != TopAbs_ON );
4278 bool ElementsOnShape::TClassifier::isOutOfBox (const gp_Pnt& p)
4280 return myBox.IsOut( p.XYZ() );
4283 bool ElementsOnShape::TClassifier::isOutOfFace (const gp_Pnt& p)
4285 myProjFace.Perform( p );
4286 if ( myProjFace.IsDone() && myProjFace.LowerDistance() <= myTol )
4288 // check relatively to the face
4289 Quantity_Parameter u, v;
4290 myProjFace.LowerDistanceParameters(u, v);
4291 gp_Pnt2d aProjPnt (u, v);
4292 BRepClass_FaceClassifier aClsf ( TopoDS::Face( myShape ), aProjPnt, myTol );
4293 if ( aClsf.State() == TopAbs_IN || aClsf.State() == TopAbs_ON )
4299 bool ElementsOnShape::TClassifier::isOutOfEdge (const gp_Pnt& p)
4301 myProjEdge.Perform( p );
4302 return ! ( myProjEdge.NbPoints() > 0 && myProjEdge.LowerDistance() <= myTol );
4305 bool ElementsOnShape::TClassifier::isOutOfVertex(const gp_Pnt& p)
4307 return ( myVertexXYZ.Distance( p ) > myTol );
4310 bool ElementsOnShape::TClassifier::isBox (const TopoDS_Shape& theShape)
4312 TopTools_IndexedMapOfShape vMap;
4313 TopExp::MapShapes( theShape, TopAbs_VERTEX, vMap );
4314 if ( vMap.Extent() != 8 )
4318 for ( int i = 1; i <= 8; ++i )
4319 myBox.Add( BRep_Tool::Pnt( TopoDS::Vertex( vMap( i ))).XYZ() );
4321 gp_XYZ pMin = myBox.CornerMin(), pMax = myBox.CornerMax();
4322 for ( int i = 1; i <= 8; ++i )
4324 gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex( vMap( i )));
4325 for ( int iC = 1; iC <= 3; ++ iC )
4327 double d1 = Abs( pMin.Coord( iC ) - p.Coord( iC ));
4328 double d2 = Abs( pMax.Coord( iC ) - p.Coord( iC ));
4329 if ( Min( d1, d2 ) > myTol )
4333 myBox.Enlarge( myTol );
4339 Class : BelongToGeom
4340 Description : Predicate for verifying whether entity belongs to
4341 specified geometrical support
4344 BelongToGeom::BelongToGeom()
4346 myType(SMDSAbs_All),
4347 myIsSubshape(false),
4348 myTolerance(Precision::Confusion())
4351 void BelongToGeom::SetMesh( const SMDS_Mesh* theMesh )
4353 myMeshDS = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
4357 void BelongToGeom::SetGeom( const TopoDS_Shape& theShape )
4363 static bool IsSubShape (const TopTools_IndexedMapOfShape& theMap,
4364 const TopoDS_Shape& theShape)
4366 if (theMap.Contains(theShape)) return true;
4368 if (theShape.ShapeType() == TopAbs_COMPOUND ||
4369 theShape.ShapeType() == TopAbs_COMPSOLID)
4371 TopoDS_Iterator anIt (theShape, Standard_True, Standard_True);
4372 for (; anIt.More(); anIt.Next())
4374 if (!IsSubShape(theMap, anIt.Value())) {
4384 void BelongToGeom::init()
4386 if (!myMeshDS || myShape.IsNull()) return;
4388 // is sub-shape of main shape?
4389 TopoDS_Shape aMainShape = myMeshDS->ShapeToMesh();
4390 if (aMainShape.IsNull()) {
4391 myIsSubshape = false;
4394 TopTools_IndexedMapOfShape aMap;
4395 TopExp::MapShapes(aMainShape, aMap);
4396 myIsSubshape = IsSubShape(aMap, myShape);
4399 //if (!myIsSubshape) // to be always ready to check an element not bound to geometry
4401 myElementsOnShapePtr.reset(new ElementsOnShape());
4402 myElementsOnShapePtr->SetTolerance(myTolerance);
4403 myElementsOnShapePtr->SetAllNodes(true); // "belong", while false means "lays on"
4404 myElementsOnShapePtr->SetMesh(myMeshDS);
4405 myElementsOnShapePtr->SetShape(myShape, myType);
4409 static bool IsContains( const SMESHDS_Mesh* theMeshDS,
4410 const TopoDS_Shape& theShape,
4411 const SMDS_MeshElement* theElem,
4412 TopAbs_ShapeEnum theFindShapeEnum,
4413 TopAbs_ShapeEnum theAvoidShapeEnum = TopAbs_SHAPE )
4415 TopExp_Explorer anExp( theShape,theFindShapeEnum,theAvoidShapeEnum );
4417 while( anExp.More() )
4419 const TopoDS_Shape& aShape = anExp.Current();
4420 if( SMESHDS_SubMesh* aSubMesh = theMeshDS->MeshElements( aShape ) ){
4421 if( aSubMesh->Contains( theElem ) )
4429 bool BelongToGeom::IsSatisfy (long theId)
4431 if (myMeshDS == 0 || myShape.IsNull())
4436 return myElementsOnShapePtr->IsSatisfy(theId);
4440 if (myType == SMDSAbs_Node)
4442 if( const SMDS_MeshNode* aNode = myMeshDS->FindNode( theId ) )
4444 if ( aNode->getshapeId() < 1 )
4445 return myElementsOnShapePtr->IsSatisfy(theId);
4447 const SMDS_PositionPtr& aPosition = aNode->GetPosition();
4448 SMDS_TypeOfPosition aTypeOfPosition = aPosition->GetTypeOfPosition();
4449 switch( aTypeOfPosition )
4451 case SMDS_TOP_VERTEX : return ( IsContains( myMeshDS,myShape,aNode,TopAbs_VERTEX ));
4452 case SMDS_TOP_EDGE : return ( IsContains( myMeshDS,myShape,aNode,TopAbs_EDGE ));
4453 case SMDS_TOP_FACE : return ( IsContains( myMeshDS,myShape,aNode,TopAbs_FACE ));
4454 case SMDS_TOP_3DSPACE: return ( IsContains( myMeshDS,myShape,aNode,TopAbs_SOLID ) ||
4455 IsContains( myMeshDS,myShape,aNode,TopAbs_SHELL ));
4461 if ( const SMDS_MeshElement* anElem = myMeshDS->FindElement( theId ))
4463 if ( anElem->getshapeId() < 1 )
4464 return myElementsOnShapePtr->IsSatisfy(theId);
4466 if( myType == SMDSAbs_All )
4468 return ( IsContains( myMeshDS,myShape,anElem,TopAbs_EDGE ) ||
4469 IsContains( myMeshDS,myShape,anElem,TopAbs_FACE ) ||
4470 IsContains( myMeshDS,myShape,anElem,TopAbs_SOLID )||
4471 IsContains( myMeshDS,myShape,anElem,TopAbs_SHELL ));
4473 else if( myType == anElem->GetType() )
4477 case SMDSAbs_Edge : return ( IsContains( myMeshDS,myShape,anElem,TopAbs_EDGE ));
4478 case SMDSAbs_Face : return ( IsContains( myMeshDS,myShape,anElem,TopAbs_FACE ));
4479 case SMDSAbs_Volume: return ( IsContains( myMeshDS,myShape,anElem,TopAbs_SOLID )||
4480 IsContains( myMeshDS,myShape,anElem,TopAbs_SHELL ));
4489 void BelongToGeom::SetType (SMDSAbs_ElementType theType)
4495 SMDSAbs_ElementType BelongToGeom::GetType() const
4500 TopoDS_Shape BelongToGeom::GetShape()
4505 const SMESHDS_Mesh* BelongToGeom::GetMeshDS() const
4510 void BelongToGeom::SetTolerance (double theTolerance)
4512 myTolerance = theTolerance;
4517 double BelongToGeom::GetTolerance()
4524 Description : Predicate for verifying whether entiy lying or partially lying on
4525 specified geometrical support
4528 LyingOnGeom::LyingOnGeom()
4530 myType(SMDSAbs_All),
4531 myIsSubshape(false),
4532 myTolerance(Precision::Confusion())
4535 void LyingOnGeom::SetMesh( const SMDS_Mesh* theMesh )
4537 myMeshDS = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
4541 void LyingOnGeom::SetGeom( const TopoDS_Shape& theShape )
4547 void LyingOnGeom::init()
4549 if (!myMeshDS || myShape.IsNull()) return;
4551 // is sub-shape of main shape?
4552 TopoDS_Shape aMainShape = myMeshDS->ShapeToMesh();
4553 if (aMainShape.IsNull()) {
4554 myIsSubshape = false;
4557 myIsSubshape = myMeshDS->IsGroupOfSubShapes( myShape );
4562 TopTools_IndexedMapOfShape shapes;
4563 TopExp::MapShapes( myShape, shapes );
4564 mySubShapesIDs.Clear();
4565 for ( int i = 1; i <= shapes.Extent(); ++i )
4567 int subID = myMeshDS->ShapeToIndex( shapes( i ));
4569 mySubShapesIDs.Add( subID );
4574 myElementsOnShapePtr.reset(new ElementsOnShape());
4575 myElementsOnShapePtr->SetTolerance(myTolerance);
4576 myElementsOnShapePtr->SetAllNodes(false); // lays on, while true means "belong"
4577 myElementsOnShapePtr->SetMesh(myMeshDS);
4578 myElementsOnShapePtr->SetShape(myShape, myType);
4582 bool LyingOnGeom::IsSatisfy( long theId )
4584 if ( myMeshDS == 0 || myShape.IsNull() )
4589 return myElementsOnShapePtr->IsSatisfy(theId);
4594 const SMDS_MeshElement* elem =
4595 ( myType == SMDSAbs_Node ) ? myMeshDS->FindNode( theId ) : myMeshDS->FindElement( theId );
4597 if ( mySubShapesIDs.Contains( elem->getshapeId() ))
4600 if ( elem->GetType() != SMDSAbs_Node )
4602 SMDS_ElemIteratorPtr nodeItr = elem->nodesIterator();
4603 while ( nodeItr->more() )
4605 const SMDS_MeshElement* aNode = nodeItr->next();
4606 if ( mySubShapesIDs.Contains( aNode->getshapeId() ))
4614 void LyingOnGeom::SetType( SMDSAbs_ElementType theType )
4620 SMDSAbs_ElementType LyingOnGeom::GetType() const
4625 TopoDS_Shape LyingOnGeom::GetShape()
4630 const SMESHDS_Mesh* LyingOnGeom::GetMeshDS() const
4635 void LyingOnGeom::SetTolerance (double theTolerance)
4637 myTolerance = theTolerance;
4642 double LyingOnGeom::GetTolerance()
4647 bool LyingOnGeom::Contains( const SMESHDS_Mesh* theMeshDS,
4648 const TopoDS_Shape& theShape,
4649 const SMDS_MeshElement* theElem,
4650 TopAbs_ShapeEnum theFindShapeEnum,
4651 TopAbs_ShapeEnum theAvoidShapeEnum )
4653 // if (IsContains(theMeshDS, theShape, theElem, theFindShapeEnum, theAvoidShapeEnum))
4656 // TopTools_MapOfShape aSubShapes;
4657 // TopExp_Explorer exp( theShape, theFindShapeEnum, theAvoidShapeEnum );
4658 // for ( ; exp.More(); exp.Next() )
4660 // const TopoDS_Shape& aShape = exp.Current();
4661 // if ( !aSubShapes.Add( aShape )) continue;
4663 // if ( SMESHDS_SubMesh* aSubMesh = theMeshDS->MeshElements( aShape ))
4665 // if ( aSubMesh->Contains( theElem ))
4668 // SMDS_ElemIteratorPtr nodeItr = theElem->nodesIterator();
4669 // while ( nodeItr->more() )
4671 // const SMDS_MeshElement* aNode = nodeItr->next();
4672 // if ( aSubMesh->Contains( aNode ))
4680 TSequenceOfXYZ::TSequenceOfXYZ(): myElem(0)
4683 TSequenceOfXYZ::TSequenceOfXYZ(size_type n) : myArray(n), myElem(0)
4686 TSequenceOfXYZ::TSequenceOfXYZ(size_type n, const gp_XYZ& t) : myArray(n,t), myElem(0)
4689 TSequenceOfXYZ::TSequenceOfXYZ(const TSequenceOfXYZ& theSequenceOfXYZ) : myArray(theSequenceOfXYZ.myArray), myElem(theSequenceOfXYZ.myElem)
4692 template <class InputIterator>
4693 TSequenceOfXYZ::TSequenceOfXYZ(InputIterator theBegin, InputIterator theEnd): myArray(theBegin,theEnd), myElem(0)
4696 TSequenceOfXYZ::~TSequenceOfXYZ()
4699 TSequenceOfXYZ& TSequenceOfXYZ::operator=(const TSequenceOfXYZ& theSequenceOfXYZ)
4701 myArray = theSequenceOfXYZ.myArray;
4702 myElem = theSequenceOfXYZ.myElem;
4706 gp_XYZ& TSequenceOfXYZ::operator()(size_type n)
4708 return myArray[n-1];
4711 const gp_XYZ& TSequenceOfXYZ::operator()(size_type n) const
4713 return myArray[n-1];
4716 void TSequenceOfXYZ::clear()
4721 void TSequenceOfXYZ::reserve(size_type n)
4726 void TSequenceOfXYZ::push_back(const gp_XYZ& v)
4728 myArray.push_back(v);
4731 TSequenceOfXYZ::size_type TSequenceOfXYZ::size() const
4733 return myArray.size();
4736 SMDSAbs_EntityType TSequenceOfXYZ::getElementEntity() const
4738 return myElem ? myElem->GetEntityType() : SMDSEntity_Last;
4741 TMeshModifTracer::TMeshModifTracer():
4742 myMeshModifTime(0), myMesh(0)
4745 void TMeshModifTracer::SetMesh( const SMDS_Mesh* theMesh )
4747 if ( theMesh != myMesh )
4748 myMeshModifTime = 0;
4751 bool TMeshModifTracer::IsMeshModified()
4753 bool modified = false;
4756 modified = ( myMeshModifTime != myMesh->GetMTime() );
4757 myMeshModifTime = myMesh->GetMTime();