1 // Copyright (C) 2007-2012 CEA/DEN, EDF R&D, OPEN CASCADE
3 // Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License.
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
23 #include "SMESH_ControlsDef.hxx"
28 #include <BRepAdaptor_Surface.hxx>
29 #include <BRepClass_FaceClassifier.hxx>
30 #include <BRep_Tool.hxx>
34 #include <TopoDS_Edge.hxx>
35 #include <TopoDS_Face.hxx>
36 #include <TopoDS_Shape.hxx>
37 #include <TopoDS_Vertex.hxx>
38 #include <TopoDS_Iterator.hxx>
40 #include <Geom_CylindricalSurface.hxx>
41 #include <Geom_Plane.hxx>
42 #include <Geom_Surface.hxx>
44 #include <Precision.hxx>
45 #include <TColStd_MapIteratorOfMapOfInteger.hxx>
46 #include <TColStd_MapOfInteger.hxx>
47 #include <TColStd_SequenceOfAsciiString.hxx>
48 #include <TColgp_Array1OfXYZ.hxx>
51 #include <gp_Cylinder.hxx>
58 #include "SMDS_Mesh.hxx"
59 #include "SMDS_Iterator.hxx"
60 #include "SMDS_MeshElement.hxx"
61 #include "SMDS_MeshNode.hxx"
62 #include "SMDS_VolumeTool.hxx"
63 #include "SMDS_QuadraticFaceOfNodes.hxx"
64 #include "SMDS_QuadraticEdge.hxx"
66 #include "SMESHDS_Mesh.hxx"
67 #include "SMESHDS_GroupBase.hxx"
69 #include "SMESH_OctreeNode.hxx"
71 #include <vtkMeshQuality.h>
79 inline gp_XYZ gpXYZ(const SMDS_MeshNode* aNode )
81 return gp_XYZ(aNode->X(), aNode->Y(), aNode->Z() );
84 inline double getAngle( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
86 gp_Vec v1( P1 - P2 ), v2( P3 - P2 );
88 return v1.Magnitude() < gp::Resolution() ||
89 v2.Magnitude() < gp::Resolution() ? 0 : v1.Angle( v2 );
92 inline double getArea( const gp_XYZ& P1, const gp_XYZ& P2, const gp_XYZ& P3 )
94 gp_Vec aVec1( P2 - P1 );
95 gp_Vec aVec2( P3 - P1 );
96 return ( aVec1 ^ aVec2 ).Magnitude() * 0.5;
99 inline double getArea( const gp_Pnt& P1, const gp_Pnt& P2, const gp_Pnt& P3 )
101 return getArea( P1.XYZ(), P2.XYZ(), P3.XYZ() );
106 inline double getDistance( const gp_XYZ& P1, const gp_XYZ& P2 )
108 double aDist = gp_Pnt( P1 ).Distance( gp_Pnt( P2 ) );
112 int getNbMultiConnection( const SMDS_Mesh* theMesh, const int theId )
117 const SMDS_MeshElement* anEdge = theMesh->FindElement( theId );
118 if ( anEdge == 0 || anEdge->GetType() != SMDSAbs_Edge/* || anEdge->NbNodes() != 2 */)
121 // for each pair of nodes in anEdge (there are 2 pairs in a quadratic edge)
122 // count elements containing both nodes of the pair.
123 // Note that there may be such cases for a quadratic edge (a horizontal line):
128 // +-----+------+ +-----+------+
131 // result sould be 2 in both cases
133 int aResult0 = 0, aResult1 = 0;
134 // last node, it is a medium one in a quadratic edge
135 const SMDS_MeshNode* aLastNode = anEdge->GetNode( anEdge->NbNodes() - 1 );
136 const SMDS_MeshNode* aNode0 = anEdge->GetNode( 0 );
137 const SMDS_MeshNode* aNode1 = anEdge->GetNode( 1 );
138 if ( aNode1 == aLastNode ) aNode1 = 0;
140 SMDS_ElemIteratorPtr anElemIter = aLastNode->GetInverseElementIterator();
141 while( anElemIter->more() ) {
142 const SMDS_MeshElement* anElem = anElemIter->next();
143 if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
144 SMDS_ElemIteratorPtr anIter = anElem->nodesIterator();
145 while ( anIter->more() ) {
146 if ( const SMDS_MeshElement* anElemNode = anIter->next() ) {
147 if ( anElemNode == aNode0 ) {
149 if ( !aNode1 ) break; // not a quadratic edge
151 else if ( anElemNode == aNode1 )
157 int aResult = std::max ( aResult0, aResult1 );
159 // TColStd_MapOfInteger aMap;
161 // SMDS_ElemIteratorPtr anIter = anEdge->nodesIterator();
162 // if ( anIter != 0 ) {
163 // while( anIter->more() ) {
164 // const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
167 // SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
168 // while( anElemIter->more() ) {
169 // const SMDS_MeshElement* anElem = anElemIter->next();
170 // if ( anElem != 0 && anElem->GetType() != SMDSAbs_Edge ) {
171 // int anId = anElem->GetID();
173 // if ( anIter->more() ) // i.e. first node
175 // else if ( aMap.Contains( anId ) )
185 gp_XYZ getNormale( const SMDS_MeshFace* theFace, bool* ok=0 )
187 int aNbNode = theFace->NbNodes();
189 gp_XYZ q1 = gpXYZ( theFace->GetNode(1)) - gpXYZ( theFace->GetNode(0));
190 gp_XYZ q2 = gpXYZ( theFace->GetNode(2)) - gpXYZ( theFace->GetNode(0));
193 gp_XYZ q3 = gpXYZ( theFace->GetNode(3)) - gpXYZ( theFace->GetNode(0));
196 double len = n.Modulus();
197 bool zeroLen = ( len <= numeric_limits<double>::min());
201 if (ok) *ok = !zeroLen;
209 using namespace SMESH::Controls;
216 Class : NumericalFunctor
217 Description : Base class for numerical functors
219 NumericalFunctor::NumericalFunctor():
225 void NumericalFunctor::SetMesh( const SMDS_Mesh* theMesh )
230 bool NumericalFunctor::GetPoints(const int theId,
231 TSequenceOfXYZ& theRes ) const
238 return GetPoints( myMesh->FindElement( theId ), theRes );
241 bool NumericalFunctor::GetPoints(const SMDS_MeshElement* anElem,
242 TSequenceOfXYZ& theRes )
249 theRes.reserve( anElem->NbNodes() );
251 // Get nodes of the element
252 SMDS_ElemIteratorPtr anIter;
254 if ( anElem->IsQuadratic() ) {
255 switch ( anElem->GetType() ) {
257 anIter = dynamic_cast<const SMDS_VtkEdge*>
258 (anElem)->interlacedNodesElemIterator();
261 anIter = dynamic_cast<const SMDS_VtkFace*>
262 (anElem)->interlacedNodesElemIterator();
265 anIter = anElem->nodesIterator();
270 anIter = anElem->nodesIterator();
274 while( anIter->more() ) {
275 if ( const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>( anIter->next() ))
276 theRes.push_back( gp_XYZ( aNode->X(), aNode->Y(), aNode->Z() ) );
283 long NumericalFunctor::GetPrecision() const
288 void NumericalFunctor::SetPrecision( const long thePrecision )
290 myPrecision = thePrecision;
291 myPrecisionValue = pow( 10., (double)( myPrecision ) );
294 double NumericalFunctor::GetValue( long theId )
298 myCurrElement = myMesh->FindElement( theId );
301 if ( GetPoints( theId, P ))
302 aVal = Round( GetValue( P ));
307 double NumericalFunctor::Round( const double & aVal )
309 return ( myPrecision >= 0 ) ? floor( aVal * myPrecisionValue + 0.5 ) / myPrecisionValue : aVal;
312 //================================================================================
314 * \brief Return histogram of functor values
315 * \param nbIntervals - number of intervals
316 * \param nbEvents - number of mesh elements having values within i-th interval
317 * \param funValues - boundaries of intervals
318 * \param elements - elements to check vulue of; empty list means "of all"
319 * \param minmax - boundaries of diapason of values to divide into intervals
321 //================================================================================
323 void NumericalFunctor::GetHistogram(int nbIntervals,
324 std::vector<int>& nbEvents,
325 std::vector<double>& funValues,
326 const vector<int>& elements,
327 const double* minmax)
329 if ( nbIntervals < 1 ||
331 !myMesh->GetMeshInfo().NbElements( GetType() ))
333 nbEvents.resize( nbIntervals, 0 );
334 funValues.resize( nbIntervals+1 );
336 // get all values sorted
337 std::multiset< double > values;
338 if ( elements.empty() )
340 SMDS_ElemIteratorPtr elemIt = myMesh->elementsIterator(GetType());
341 while ( elemIt->more() )
342 values.insert( GetValue( elemIt->next()->GetID() ));
346 vector<int>::const_iterator id = elements.begin();
347 for ( ; id != elements.end(); ++id )
348 values.insert( GetValue( *id ));
353 funValues[0] = minmax[0];
354 funValues[nbIntervals] = minmax[1];
358 funValues[0] = *values.begin();
359 funValues[nbIntervals] = *values.rbegin();
361 // case nbIntervals == 1
362 if ( nbIntervals == 1 )
364 nbEvents[0] = values.size();
368 if (funValues.front() == funValues.back())
370 nbEvents.resize( 1 );
371 nbEvents[0] = values.size();
372 funValues[1] = funValues.back();
373 funValues.resize( 2 );
376 std::multiset< double >::iterator min = values.begin(), max;
377 for ( int i = 0; i < nbIntervals; ++i )
379 // find end value of i-th interval
380 double r = (i+1) / double( nbIntervals );
381 funValues[i+1] = funValues.front() * (1-r) + funValues.back() * r;
383 // count values in the i-th interval if there are any
384 if ( min != values.end() && *min <= funValues[i+1] )
386 // find the first value out of the interval
387 max = values.upper_bound( funValues[i+1] ); // max is greater than funValues[i+1], or end()
388 nbEvents[i] = std::distance( min, max );
392 // add values larger than minmax[1]
393 nbEvents.back() += std::distance( min, values.end() );
396 //=======================================================================
397 //function : GetValue
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 //=======================================================================
412 //function : GetBadRate
413 //purpose : meaningless as it is not quality control functor
414 //=======================================================================
416 double Volume::GetBadRate( double Value, int /*nbNodes*/ ) const
421 //=======================================================================
424 //=======================================================================
426 SMDSAbs_ElementType Volume::GetType() const
428 return SMDSAbs_Volume;
433 Class : MaxElementLength2D
434 Description : Functor calculating maximum length of 2D element
437 double MaxElementLength2D::GetValue( long theElementId )
440 if( GetPoints( theElementId, P ) ) {
442 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
443 SMDSAbs_ElementType aType = aElem->GetType();
447 if( len == 3 ) { // triangles
448 double L1 = getDistance(P( 1 ),P( 2 ));
449 double L2 = getDistance(P( 2 ),P( 3 ));
450 double L3 = getDistance(P( 3 ),P( 1 ));
451 aVal = Max(L1,Max(L2,L3));
454 else if( len == 4 ) { // quadrangles
455 double L1 = getDistance(P( 1 ),P( 2 ));
456 double L2 = getDistance(P( 2 ),P( 3 ));
457 double L3 = getDistance(P( 3 ),P( 4 ));
458 double L4 = getDistance(P( 4 ),P( 1 ));
459 double D1 = getDistance(P( 1 ),P( 3 ));
460 double D2 = getDistance(P( 2 ),P( 4 ));
461 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
464 else if( len == 6 ) { // quadratic triangles
465 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
466 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
467 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
468 aVal = Max(L1,Max(L2,L3));
471 else if( len == 8 || len == 9 ) { // quadratic quadrangles
472 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
473 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
474 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
475 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
476 double D1 = getDistance(P( 1 ),P( 5 ));
477 double D2 = getDistance(P( 3 ),P( 7 ));
478 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(D1,D2));
483 if( myPrecision >= 0 )
485 double prec = pow( 10., (double)myPrecision );
486 aVal = floor( aVal * prec + 0.5 ) / prec;
493 double MaxElementLength2D::GetBadRate( double Value, int /*nbNodes*/ ) const
498 SMDSAbs_ElementType MaxElementLength2D::GetType() const
504 Class : MaxElementLength3D
505 Description : Functor calculating maximum length of 3D element
508 double MaxElementLength3D::GetValue( long theElementId )
511 if( GetPoints( theElementId, P ) ) {
513 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
514 SMDSAbs_ElementType aType = aElem->GetType();
518 if( len == 4 ) { // tetras
519 double L1 = getDistance(P( 1 ),P( 2 ));
520 double L2 = getDistance(P( 2 ),P( 3 ));
521 double L3 = getDistance(P( 3 ),P( 1 ));
522 double L4 = getDistance(P( 1 ),P( 4 ));
523 double L5 = getDistance(P( 2 ),P( 4 ));
524 double L6 = getDistance(P( 3 ),P( 4 ));
525 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
528 else if( len == 5 ) { // pyramids
529 double L1 = getDistance(P( 1 ),P( 2 ));
530 double L2 = getDistance(P( 2 ),P( 3 ));
531 double L3 = getDistance(P( 3 ),P( 4 ));
532 double L4 = getDistance(P( 4 ),P( 1 ));
533 double L5 = getDistance(P( 1 ),P( 5 ));
534 double L6 = getDistance(P( 2 ),P( 5 ));
535 double L7 = getDistance(P( 3 ),P( 5 ));
536 double L8 = getDistance(P( 4 ),P( 5 ));
537 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
538 aVal = Max(aVal,Max(L7,L8));
541 else if( len == 6 ) { // pentas
542 double L1 = getDistance(P( 1 ),P( 2 ));
543 double L2 = getDistance(P( 2 ),P( 3 ));
544 double L3 = getDistance(P( 3 ),P( 1 ));
545 double L4 = getDistance(P( 4 ),P( 5 ));
546 double L5 = getDistance(P( 5 ),P( 6 ));
547 double L6 = getDistance(P( 6 ),P( 4 ));
548 double L7 = getDistance(P( 1 ),P( 4 ));
549 double L8 = getDistance(P( 2 ),P( 5 ));
550 double L9 = getDistance(P( 3 ),P( 6 ));
551 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
552 aVal = Max(aVal,Max(Max(L7,L8),L9));
555 else if( len == 8 ) { // hexas
556 double L1 = getDistance(P( 1 ),P( 2 ));
557 double L2 = getDistance(P( 2 ),P( 3 ));
558 double L3 = getDistance(P( 3 ),P( 4 ));
559 double L4 = getDistance(P( 4 ),P( 1 ));
560 double L5 = getDistance(P( 5 ),P( 6 ));
561 double L6 = getDistance(P( 6 ),P( 7 ));
562 double L7 = getDistance(P( 7 ),P( 8 ));
563 double L8 = getDistance(P( 8 ),P( 5 ));
564 double L9 = getDistance(P( 1 ),P( 5 ));
565 double L10= getDistance(P( 2 ),P( 6 ));
566 double L11= getDistance(P( 3 ),P( 7 ));
567 double L12= getDistance(P( 4 ),P( 8 ));
568 double D1 = getDistance(P( 1 ),P( 7 ));
569 double D2 = getDistance(P( 2 ),P( 8 ));
570 double D3 = getDistance(P( 3 ),P( 5 ));
571 double D4 = getDistance(P( 4 ),P( 6 ));
572 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
573 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
574 aVal = Max(aVal,Max(L11,L12));
575 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
578 else if( len == 12 ) { // hexagonal prism
579 for ( int i1 = 1; i1 < 12; ++i1 )
580 for ( int i2 = i1+1; i1 <= 12; ++i1 )
581 aVal = Max( aVal, getDistance(P( i1 ),P( i2 )));
584 else if( len == 10 ) { // quadratic tetras
585 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
586 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
587 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
588 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
589 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
590 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
591 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
594 else if( len == 13 ) { // quadratic pyramids
595 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
596 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
597 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
598 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
599 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
600 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
601 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
602 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
603 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
604 aVal = Max(aVal,Max(L7,L8));
607 else if( len == 15 ) { // quadratic pentas
608 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
609 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
610 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
611 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
612 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
613 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
614 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
615 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
616 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
617 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
618 aVal = Max(aVal,Max(Max(L7,L8),L9));
621 else if( len == 20 || len == 27 ) { // quadratic hexas
622 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
623 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
624 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
625 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
626 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
627 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
628 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
629 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
630 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
631 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
632 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
633 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
634 double D1 = getDistance(P( 1 ),P( 7 ));
635 double D2 = getDistance(P( 2 ),P( 8 ));
636 double D3 = getDistance(P( 3 ),P( 5 ));
637 double D4 = getDistance(P( 4 ),P( 6 ));
638 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
639 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
640 aVal = Max(aVal,Max(L11,L12));
641 aVal = Max(aVal,Max(Max(D1,D2),Max(D3,D4)));
644 else if( len > 1 && aElem->IsPoly() ) { // polys
645 // get the maximum distance between all pairs of nodes
646 for( int i = 1; i <= len; i++ ) {
647 for( int j = 1; j <= len; j++ ) {
648 if( j > i ) { // optimization of the loop
649 double D = getDistance( P(i), P(j) );
650 aVal = Max( aVal, D );
657 if( myPrecision >= 0 )
659 double prec = pow( 10., (double)myPrecision );
660 aVal = floor( aVal * prec + 0.5 ) / prec;
667 double MaxElementLength3D::GetBadRate( double Value, int /*nbNodes*/ ) const
672 SMDSAbs_ElementType MaxElementLength3D::GetType() const
674 return SMDSAbs_Volume;
680 Description : Functor for calculation of minimum angle
683 double MinimumAngle::GetValue( const TSequenceOfXYZ& P )
690 aMin = getAngle(P( P.size() ), P( 1 ), P( 2 ));
691 aMin = Min(aMin,getAngle(P( P.size()-1 ), P( P.size() ), P( 1 )));
693 for (int i=2; i<P.size();i++){
694 double A0 = getAngle( P( i-1 ), P( i ), P( i+1 ) );
698 return aMin * 180.0 / M_PI;
701 double MinimumAngle::GetBadRate( double Value, int nbNodes ) const
703 //const double aBestAngle = PI / nbNodes;
704 const double aBestAngle = 180.0 - ( 360.0 / double(nbNodes) );
705 return ( fabs( aBestAngle - Value ));
708 SMDSAbs_ElementType MinimumAngle::GetType() const
716 Description : Functor for calculating aspect ratio
718 double AspectRatio::GetValue( long theId )
721 myCurrElement = myMesh->FindElement( theId );
722 if ( myCurrElement && myCurrElement->GetVtkType() == VTK_QUAD )
725 vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myCurrElement->getMeshId()]->getGrid();
726 if ( vtkCell* avtkCell = grid->GetCell( myCurrElement->getVtkId() ))
727 aVal = Round( vtkMeshQuality::QuadAspectRatio( avtkCell ));
732 if ( GetPoints( myCurrElement, P ))
733 aVal = Round( GetValue( P ));
738 double AspectRatio::GetValue( const TSequenceOfXYZ& P )
740 // According to "Mesh quality control" by Nadir Bouhamau referring to
741 // Pascal Jean Frey and Paul-Louis George. Maillages, applications aux elements finis.
742 // Hermes Science publications, Paris 1999 ISBN 2-7462-0024-4
745 int nbNodes = P.size();
750 // Compute aspect ratio
752 if ( nbNodes == 3 ) {
753 // Compute lengths of the sides
754 std::vector< double > aLen (nbNodes);
755 for ( int i = 0; i < nbNodes - 1; i++ )
756 aLen[ i ] = getDistance( P( i + 1 ), P( i + 2 ) );
757 aLen[ nbNodes - 1 ] = getDistance( P( 1 ), P( nbNodes ) );
758 // Q = alfa * h * p / S, where
760 // alfa = sqrt( 3 ) / 6
761 // h - length of the longest edge
762 // p - half perimeter
763 // S - triangle surface
764 const double alfa = sqrt( 3. ) / 6.;
765 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
766 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
767 double anArea = getArea( P( 1 ), P( 2 ), P( 3 ) );
768 if ( anArea <= Precision::Confusion() )
770 return alfa * maxLen * half_perimeter / anArea;
772 else if ( nbNodes == 6 ) { // quadratic triangles
773 // Compute lengths of the sides
774 std::vector< double > aLen (3);
775 aLen[0] = getDistance( P(1), P(3) );
776 aLen[1] = getDistance( P(3), P(5) );
777 aLen[2] = getDistance( P(5), P(1) );
778 // Q = alfa * h * p / S, where
780 // alfa = sqrt( 3 ) / 6
781 // h - length of the longest edge
782 // p - half perimeter
783 // S - triangle surface
784 const double alfa = sqrt( 3. ) / 6.;
785 double maxLen = Max( aLen[ 0 ], Max( aLen[ 1 ], aLen[ 2 ] ) );
786 double half_perimeter = ( aLen[0] + aLen[1] + aLen[2] ) / 2.;
787 double anArea = getArea( P(1), P(3), P(5) );
788 if ( anArea <= Precision::Confusion() )
790 return alfa * maxLen * half_perimeter / anArea;
792 else if( nbNodes == 4 ) { // quadrangle
793 // Compute lengths of the sides
794 std::vector< double > aLen (4);
795 aLen[0] = getDistance( P(1), P(2) );
796 aLen[1] = getDistance( P(2), P(3) );
797 aLen[2] = getDistance( P(3), P(4) );
798 aLen[3] = getDistance( P(4), P(1) );
799 // Compute lengths of the diagonals
800 std::vector< double > aDia (2);
801 aDia[0] = getDistance( P(1), P(3) );
802 aDia[1] = getDistance( P(2), P(4) );
803 // Compute areas of all triangles which can be built
804 // taking three nodes of the quadrangle
805 std::vector< double > anArea (4);
806 anArea[0] = getArea( P(1), P(2), P(3) );
807 anArea[1] = getArea( P(1), P(2), P(4) );
808 anArea[2] = getArea( P(1), P(3), P(4) );
809 anArea[3] = getArea( P(2), P(3), P(4) );
810 // Q = alpha * L * C1 / C2, where
812 // alpha = sqrt( 1/32 )
813 // L = max( L1, L2, L3, L4, D1, D2 )
814 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
815 // C2 = min( S1, S2, S3, S4 )
816 // Li - lengths of the edges
817 // Di - lengths of the diagonals
818 // Si - areas of the triangles
819 const double alpha = sqrt( 1 / 32. );
820 double L = Max( aLen[ 0 ],
824 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
825 double C1 = sqrt( ( aLen[0] * aLen[0] +
828 aLen[3] * aLen[3] ) / 4. );
829 double C2 = Min( anArea[ 0 ],
831 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
832 if ( C2 <= Precision::Confusion() )
834 return alpha * L * C1 / C2;
836 else if( nbNodes == 8 || nbNodes == 9 ) { // nbNodes==8 - quadratic quadrangle
837 // Compute lengths of the sides
838 std::vector< double > aLen (4);
839 aLen[0] = getDistance( P(1), P(3) );
840 aLen[1] = getDistance( P(3), P(5) );
841 aLen[2] = getDistance( P(5), P(7) );
842 aLen[3] = getDistance( P(7), P(1) );
843 // Compute lengths of the diagonals
844 std::vector< double > aDia (2);
845 aDia[0] = getDistance( P(1), P(5) );
846 aDia[1] = getDistance( P(3), P(7) );
847 // Compute areas of all triangles which can be built
848 // taking three nodes of the quadrangle
849 std::vector< double > anArea (4);
850 anArea[0] = getArea( P(1), P(3), P(5) );
851 anArea[1] = getArea( P(1), P(3), P(7) );
852 anArea[2] = getArea( P(1), P(5), P(7) );
853 anArea[3] = getArea( P(3), P(5), P(7) );
854 // Q = alpha * L * C1 / C2, where
856 // alpha = sqrt( 1/32 )
857 // L = max( L1, L2, L3, L4, D1, D2 )
858 // C1 = sqrt( ( L1^2 + L1^2 + L1^2 + L1^2 ) / 4 )
859 // C2 = min( S1, S2, S3, S4 )
860 // Li - lengths of the edges
861 // Di - lengths of the diagonals
862 // Si - areas of the triangles
863 const double alpha = sqrt( 1 / 32. );
864 double L = Max( aLen[ 0 ],
868 Max( aDia[ 0 ], aDia[ 1 ] ) ) ) ) );
869 double C1 = sqrt( ( aLen[0] * aLen[0] +
872 aLen[3] * aLen[3] ) / 4. );
873 double C2 = Min( anArea[ 0 ],
875 Min( anArea[ 2 ], anArea[ 3 ] ) ) );
876 if ( C2 <= Precision::Confusion() )
878 return alpha * L * C1 / C2;
883 double AspectRatio::GetBadRate( double Value, int /*nbNodes*/ ) const
885 // the aspect ratio is in the range [1.0,infinity]
886 // < 1.0 = very bad, zero area
889 return ( Value < 0.9 ) ? 1000 : Value / 1000.;
892 SMDSAbs_ElementType AspectRatio::GetType() const
899 Class : AspectRatio3D
900 Description : Functor for calculating aspect ratio
904 inline double getHalfPerimeter(double theTria[3]){
905 return (theTria[0] + theTria[1] + theTria[2])/2.0;
908 inline double getArea(double theHalfPerim, double theTria[3]){
909 return sqrt(theHalfPerim*
910 (theHalfPerim-theTria[0])*
911 (theHalfPerim-theTria[1])*
912 (theHalfPerim-theTria[2]));
915 inline double getVolume(double theLen[6]){
916 double a2 = theLen[0]*theLen[0];
917 double b2 = theLen[1]*theLen[1];
918 double c2 = theLen[2]*theLen[2];
919 double d2 = theLen[3]*theLen[3];
920 double e2 = theLen[4]*theLen[4];
921 double f2 = theLen[5]*theLen[5];
922 double P = 4.0*a2*b2*d2;
923 double Q = a2*(b2+d2-e2)-b2*(a2+d2-f2)-d2*(a2+b2-c2);
924 double R = (b2+d2-e2)*(a2+d2-f2)*(a2+d2-f2);
925 return sqrt(P-Q+R)/12.0;
928 inline double getVolume2(double theLen[6]){
929 double a2 = theLen[0]*theLen[0];
930 double b2 = theLen[1]*theLen[1];
931 double c2 = theLen[2]*theLen[2];
932 double d2 = theLen[3]*theLen[3];
933 double e2 = theLen[4]*theLen[4];
934 double f2 = theLen[5]*theLen[5];
936 double P = a2*e2*(b2+c2+d2+f2-a2-e2);
937 double Q = b2*f2*(a2+c2+d2+e2-b2-f2);
938 double R = c2*d2*(a2+b2+e2+f2-c2-d2);
939 double S = a2*b2*d2+b2*c2*e2+a2*c2*f2+d2*e2*f2;
941 return sqrt(P+Q+R-S)/12.0;
944 inline double getVolume(const TSequenceOfXYZ& P){
945 gp_Vec aVec1( P( 2 ) - P( 1 ) );
946 gp_Vec aVec2( P( 3 ) - P( 1 ) );
947 gp_Vec aVec3( P( 4 ) - P( 1 ) );
948 gp_Vec anAreaVec( aVec1 ^ aVec2 );
949 return fabs(aVec3 * anAreaVec) / 6.0;
952 inline double getMaxHeight(double theLen[6])
954 double aHeight = std::max(theLen[0],theLen[1]);
955 aHeight = std::max(aHeight,theLen[2]);
956 aHeight = std::max(aHeight,theLen[3]);
957 aHeight = std::max(aHeight,theLen[4]);
958 aHeight = std::max(aHeight,theLen[5]);
964 double AspectRatio3D::GetValue( long theId )
967 myCurrElement = myMesh->FindElement( theId );
968 if ( myCurrElement && myCurrElement->GetVtkType() == VTK_TETRA )
970 // Action from CoTech | ACTION 31.3:
971 // EURIWARE BO: Homogenize the formulas used to calculate the Controls in SMESH to fit with
972 // those of ParaView. The library used by ParaView for those calculations can be reused in SMESH.
973 vtkUnstructuredGrid* grid = SMDS_Mesh::_meshList[myCurrElement->getMeshId()]->getGrid();
974 if ( vtkCell* avtkCell = grid->GetCell( myCurrElement->getVtkId() ))
975 aVal = Round( vtkMeshQuality::TetAspectRatio( avtkCell ));
980 if ( GetPoints( myCurrElement, P ))
981 aVal = Round( GetValue( P ));
986 double AspectRatio3D::GetValue( const TSequenceOfXYZ& P )
988 double aQuality = 0.0;
989 if(myCurrElement->IsPoly()) return aQuality;
991 int nbNodes = P.size();
993 if(myCurrElement->IsQuadratic()) {
994 if(nbNodes==10) nbNodes=4; // quadratic tetrahedron
995 else if(nbNodes==13) nbNodes=5; // quadratic pyramid
996 else if(nbNodes==15) nbNodes=6; // quadratic pentahedron
997 else if(nbNodes==20) nbNodes=8; // quadratic hexahedron
998 else if(nbNodes==27) nbNodes=8; // quadratic hexahedron
999 else return aQuality;
1005 getDistance(P( 1 ),P( 2 )), // a
1006 getDistance(P( 2 ),P( 3 )), // b
1007 getDistance(P( 3 ),P( 1 )), // c
1008 getDistance(P( 2 ),P( 4 )), // d
1009 getDistance(P( 3 ),P( 4 )), // e
1010 getDistance(P( 1 ),P( 4 )) // f
1012 double aTria[4][3] = {
1013 {aLen[0],aLen[1],aLen[2]}, // abc
1014 {aLen[0],aLen[3],aLen[5]}, // adf
1015 {aLen[1],aLen[3],aLen[4]}, // bde
1016 {aLen[2],aLen[4],aLen[5]} // cef
1018 double aSumArea = 0.0;
1019 double aHalfPerimeter = getHalfPerimeter(aTria[0]);
1020 double anArea = getArea(aHalfPerimeter,aTria[0]);
1022 aHalfPerimeter = getHalfPerimeter(aTria[1]);
1023 anArea = getArea(aHalfPerimeter,aTria[1]);
1025 aHalfPerimeter = getHalfPerimeter(aTria[2]);
1026 anArea = getArea(aHalfPerimeter,aTria[2]);
1028 aHalfPerimeter = getHalfPerimeter(aTria[3]);
1029 anArea = getArea(aHalfPerimeter,aTria[3]);
1031 double aVolume = getVolume(P);
1032 //double aVolume = getVolume(aLen);
1033 double aHeight = getMaxHeight(aLen);
1034 static double aCoeff = sqrt(2.0)/12.0;
1035 if ( aVolume > DBL_MIN )
1036 aQuality = aCoeff*aHeight*aSumArea/aVolume;
1041 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 3 ),P( 5 )};
1042 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1045 gp_XYZ aXYZ[4] = {P( 1 ),P( 3 ),P( 4 ),P( 5 )};
1046 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1049 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 5 )};
1050 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1053 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 4 ),P( 5 )};
1054 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1060 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 6 )};
1061 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1064 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 3 )};
1065 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1068 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 6 )};
1069 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1072 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1073 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1076 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 6 )};
1077 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1080 gp_XYZ aXYZ[4] = {P( 2 ),P( 5 ),P( 4 ),P( 3 )};
1081 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1087 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 3 )};
1088 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1091 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 4 )};
1092 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1095 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 7 )};
1096 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1099 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 5 ),P( 8 )};
1100 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1103 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 3 )};
1104 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1107 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 4 )};
1108 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1111 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 7 )};
1112 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1115 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 6 ),P( 8 )};
1116 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1119 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 3 )};
1120 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1123 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 4 )};
1124 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1127 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 7 )};
1128 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1131 gp_XYZ aXYZ[4] = {P( 2 ),P( 6 ),P( 5 ),P( 8 )};
1132 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1135 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 1 )};
1136 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1139 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 2 )};
1140 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1143 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 5 )};
1144 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1147 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 8 ),P( 6 )};
1148 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1151 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 1 )};
1152 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1155 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 2 )};
1156 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1159 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 5 )};
1160 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1163 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 7 ),P( 6 )};
1164 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1167 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 1 )};
1168 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1171 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1172 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1175 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 5 )};
1176 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1179 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 6 )};
1180 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1183 gp_XYZ aXYZ[4] = {P( 4 ),P( 8 ),P( 7 ),P( 2 )};
1184 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1187 gp_XYZ aXYZ[4] = {P( 4 ),P( 5 ),P( 8 ),P( 2 )};
1188 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1191 gp_XYZ aXYZ[4] = {P( 1 ),P( 4 ),P( 5 ),P( 3 )};
1192 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1195 gp_XYZ aXYZ[4] = {P( 3 ),P( 6 ),P( 7 ),P( 1 )};
1196 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1199 gp_XYZ aXYZ[4] = {P( 2 ),P( 3 ),P( 6 ),P( 4 )};
1200 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1203 gp_XYZ aXYZ[4] = {P( 5 ),P( 6 ),P( 8 ),P( 3 )};
1204 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1207 gp_XYZ aXYZ[4] = {P( 7 ),P( 8 ),P( 6 ),P( 1 )};
1208 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1211 gp_XYZ aXYZ[4] = {P( 1 ),P( 2 ),P( 4 ),P( 7 )};
1212 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1215 gp_XYZ aXYZ[4] = {P( 3 ),P( 4 ),P( 2 ),P( 5 )};
1216 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[4])),aQuality);
1222 gp_XYZ aXYZ[8] = {P( 1 ),P( 2 ),P( 4 ),P( 5 ),P( 7 ),P( 8 ),P( 10 ),P( 11 )};
1223 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1226 gp_XYZ aXYZ[8] = {P( 2 ),P( 3 ),P( 5 ),P( 6 ),P( 8 ),P( 9 ),P( 11 ),P( 12 )};
1227 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1230 gp_XYZ aXYZ[8] = {P( 3 ),P( 4 ),P( 6 ),P( 1 ),P( 9 ),P( 10 ),P( 12 ),P( 7 )};
1231 aQuality = std::max(GetValue(TSequenceOfXYZ(&aXYZ[0],&aXYZ[8])),aQuality);
1234 } // switch(nbNodes)
1236 if ( nbNodes > 4 ) {
1237 // avaluate aspect ratio of quadranle faces
1238 AspectRatio aspect2D;
1239 SMDS_VolumeTool::VolumeType type = SMDS_VolumeTool::GetType( nbNodes );
1240 int nbFaces = SMDS_VolumeTool::NbFaces( type );
1241 TSequenceOfXYZ points(4);
1242 for ( int i = 0; i < nbFaces; ++i ) { // loop on faces of a volume
1243 if ( SMDS_VolumeTool::NbFaceNodes( type, i ) != 4 )
1245 const int* pInd = SMDS_VolumeTool::GetFaceNodesIndices( type, i, true );
1246 for ( int p = 0; p < 4; ++p ) // loop on nodes of a quadranle face
1247 points( p + 1 ) = P( pInd[ p ] + 1 );
1248 aQuality = std::max( aQuality, aspect2D.GetValue( points ));
1254 double AspectRatio3D::GetBadRate( double Value, int /*nbNodes*/ ) const
1256 // the aspect ratio is in the range [1.0,infinity]
1259 return Value / 1000.;
1262 SMDSAbs_ElementType AspectRatio3D::GetType() const
1264 return SMDSAbs_Volume;
1270 Description : Functor for calculating warping
1272 double Warping::GetValue( const TSequenceOfXYZ& P )
1274 if ( P.size() != 4 )
1277 gp_XYZ G = ( P( 1 ) + P( 2 ) + P( 3 ) + P( 4 ) ) / 4.;
1279 double A1 = ComputeA( P( 1 ), P( 2 ), P( 3 ), G );
1280 double A2 = ComputeA( P( 2 ), P( 3 ), P( 4 ), G );
1281 double A3 = ComputeA( P( 3 ), P( 4 ), P( 1 ), G );
1282 double A4 = ComputeA( P( 4 ), P( 1 ), P( 2 ), G );
1284 return Max( Max( A1, A2 ), Max( A3, A4 ) );
1287 double Warping::ComputeA( const gp_XYZ& thePnt1,
1288 const gp_XYZ& thePnt2,
1289 const gp_XYZ& thePnt3,
1290 const gp_XYZ& theG ) const
1292 double aLen1 = gp_Pnt( thePnt1 ).Distance( gp_Pnt( thePnt2 ) );
1293 double aLen2 = gp_Pnt( thePnt2 ).Distance( gp_Pnt( thePnt3 ) );
1294 double L = Min( aLen1, aLen2 ) * 0.5;
1295 if ( L < Precision::Confusion())
1298 gp_XYZ GI = ( thePnt2 + thePnt1 ) / 2. - theG;
1299 gp_XYZ GJ = ( thePnt3 + thePnt2 ) / 2. - theG;
1300 gp_XYZ N = GI.Crossed( GJ );
1302 if ( N.Modulus() < gp::Resolution() )
1307 double H = ( thePnt2 - theG ).Dot( N );
1308 return asin( fabs( H / L ) ) * 180. / M_PI;
1311 double Warping::GetBadRate( double Value, int /*nbNodes*/ ) const
1313 // the warp is in the range [0.0,PI/2]
1314 // 0.0 = good (no warp)
1315 // PI/2 = bad (face pliee)
1319 SMDSAbs_ElementType Warping::GetType() const
1321 return SMDSAbs_Face;
1327 Description : Functor for calculating taper
1329 double Taper::GetValue( const TSequenceOfXYZ& P )
1331 if ( P.size() != 4 )
1335 double J1 = getArea( P( 4 ), P( 1 ), P( 2 ) ) / 2.;
1336 double J2 = getArea( P( 3 ), P( 1 ), P( 2 ) ) / 2.;
1337 double J3 = getArea( P( 2 ), P( 3 ), P( 4 ) ) / 2.;
1338 double J4 = getArea( P( 3 ), P( 4 ), P( 1 ) ) / 2.;
1340 double JA = 0.25 * ( J1 + J2 + J3 + J4 );
1341 if ( JA <= Precision::Confusion() )
1344 double T1 = fabs( ( J1 - JA ) / JA );
1345 double T2 = fabs( ( J2 - JA ) / JA );
1346 double T3 = fabs( ( J3 - JA ) / JA );
1347 double T4 = fabs( ( J4 - JA ) / JA );
1349 return Max( Max( T1, T2 ), Max( T3, T4 ) );
1352 double Taper::GetBadRate( double Value, int /*nbNodes*/ ) const
1354 // the taper is in the range [0.0,1.0]
1355 // 0.0 = good (no taper)
1356 // 1.0 = bad (les cotes opposes sont allignes)
1360 SMDSAbs_ElementType Taper::GetType() const
1362 return SMDSAbs_Face;
1368 Description : Functor for calculating skew in degrees
1370 static inline double skewAngle( const gp_XYZ& p1, const gp_XYZ& p2, const gp_XYZ& p3 )
1372 gp_XYZ p12 = ( p2 + p1 ) / 2.;
1373 gp_XYZ p23 = ( p3 + p2 ) / 2.;
1374 gp_XYZ p31 = ( p3 + p1 ) / 2.;
1376 gp_Vec v1( p31 - p2 ), v2( p12 - p23 );
1378 return v1.Magnitude() < gp::Resolution() || v2.Magnitude() < gp::Resolution() ? 0. : v1.Angle( v2 );
1381 double Skew::GetValue( const TSequenceOfXYZ& P )
1383 if ( P.size() != 3 && P.size() != 4 )
1387 static double PI2 = M_PI / 2.;
1388 if ( P.size() == 3 )
1390 double A0 = fabs( PI2 - skewAngle( P( 3 ), P( 1 ), P( 2 ) ) );
1391 double A1 = fabs( PI2 - skewAngle( P( 1 ), P( 2 ), P( 3 ) ) );
1392 double A2 = fabs( PI2 - skewAngle( P( 2 ), P( 3 ), P( 1 ) ) );
1394 return Max( A0, Max( A1, A2 ) ) * 180. / M_PI;
1398 gp_XYZ p12 = ( P( 1 ) + P( 2 ) ) / 2.;
1399 gp_XYZ p23 = ( P( 2 ) + P( 3 ) ) / 2.;
1400 gp_XYZ p34 = ( P( 3 ) + P( 4 ) ) / 2.;
1401 gp_XYZ p41 = ( P( 4 ) + P( 1 ) ) / 2.;
1403 gp_Vec v1( p34 - p12 ), v2( p23 - p41 );
1404 double A = v1.Magnitude() <= gp::Resolution() || v2.Magnitude() <= gp::Resolution()
1405 ? 0. : fabs( PI2 - v1.Angle( v2 ) );
1408 if ( A < Precision::Angular() )
1411 return A * 180. / M_PI;
1415 double Skew::GetBadRate( double Value, int /*nbNodes*/ ) const
1417 // the skew is in the range [0.0,PI/2].
1423 SMDSAbs_ElementType Skew::GetType() const
1425 return SMDSAbs_Face;
1431 Description : Functor for calculating area
1433 double Area::GetValue( const TSequenceOfXYZ& P )
1436 if ( P.size() > 2 ) {
1437 gp_Vec aVec1( P(2) - P(1) );
1438 gp_Vec aVec2( P(3) - P(1) );
1439 gp_Vec SumVec = aVec1 ^ aVec2;
1440 for (int i=4; i<=P.size(); i++) {
1441 gp_Vec aVec1( P(i-1) - P(1) );
1442 gp_Vec aVec2( P(i) - P(1) );
1443 gp_Vec tmp = aVec1 ^ aVec2;
1446 val = SumVec.Magnitude() * 0.5;
1451 double Area::GetBadRate( double Value, int /*nbNodes*/ ) const
1453 // meaningless as it is not a quality control functor
1457 SMDSAbs_ElementType Area::GetType() const
1459 return SMDSAbs_Face;
1465 Description : Functor for calculating length of edge
1467 double Length::GetValue( const TSequenceOfXYZ& P )
1469 switch ( P.size() ) {
1470 case 2: return getDistance( P( 1 ), P( 2 ) );
1471 case 3: return getDistance( P( 1 ), P( 2 ) ) + getDistance( P( 2 ), P( 3 ) );
1476 double Length::GetBadRate( double Value, int /*nbNodes*/ ) const
1478 // meaningless as it is not quality control functor
1482 SMDSAbs_ElementType Length::GetType() const
1484 return SMDSAbs_Edge;
1489 Description : Functor for calculating length of edge
1492 double Length2D::GetValue( long theElementId)
1496 //cout<<"Length2D::GetValue"<<endl;
1497 if (GetPoints(theElementId,P)){
1498 //for(int jj=1; jj<=P.size(); jj++)
1499 // cout<<"jj="<<jj<<" P("<<P(jj).X()<<","<<P(jj).Y()<<","<<P(jj).Z()<<")"<<endl;
1501 double aVal;// = GetValue( P );
1502 const SMDS_MeshElement* aElem = myMesh->FindElement( theElementId );
1503 SMDSAbs_ElementType aType = aElem->GetType();
1512 aVal = getDistance( P( 1 ), P( 2 ) );
1515 else if (len == 3){ // quadratic edge
1516 aVal = getDistance(P( 1 ),P( 3 )) + getDistance(P( 3 ),P( 2 ));
1520 if (len == 3){ // triangles
1521 double L1 = getDistance(P( 1 ),P( 2 ));
1522 double L2 = getDistance(P( 2 ),P( 3 ));
1523 double L3 = getDistance(P( 3 ),P( 1 ));
1524 aVal = Max(L1,Max(L2,L3));
1527 else if (len == 4){ // quadrangles
1528 double L1 = getDistance(P( 1 ),P( 2 ));
1529 double L2 = getDistance(P( 2 ),P( 3 ));
1530 double L3 = getDistance(P( 3 ),P( 4 ));
1531 double L4 = getDistance(P( 4 ),P( 1 ));
1532 aVal = Max(Max(L1,L2),Max(L3,L4));
1535 if (len == 6){ // quadratic triangles
1536 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1537 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1538 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 1 ));
1539 aVal = Max(L1,Max(L2,L3));
1540 //cout<<"L1="<<L1<<" L2="<<L2<<"L3="<<L3<<" aVal="<<aVal<<endl;
1543 else if (len == 8){ // quadratic quadrangles
1544 double L1 = getDistance(P( 1 ),P( 2 )) + getDistance(P( 2 ),P( 3 ));
1545 double L2 = getDistance(P( 3 ),P( 4 )) + getDistance(P( 4 ),P( 5 ));
1546 double L3 = getDistance(P( 5 ),P( 6 )) + getDistance(P( 6 ),P( 7 ));
1547 double L4 = getDistance(P( 7 ),P( 8 )) + getDistance(P( 8 ),P( 1 ));
1548 aVal = Max(Max(L1,L2),Max(L3,L4));
1551 case SMDSAbs_Volume:
1552 if (len == 4){ // tetraidrs
1553 double L1 = getDistance(P( 1 ),P( 2 ));
1554 double L2 = getDistance(P( 2 ),P( 3 ));
1555 double L3 = getDistance(P( 3 ),P( 1 ));
1556 double L4 = getDistance(P( 1 ),P( 4 ));
1557 double L5 = getDistance(P( 2 ),P( 4 ));
1558 double L6 = getDistance(P( 3 ),P( 4 ));
1559 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1562 else if (len == 5){ // piramids
1563 double L1 = getDistance(P( 1 ),P( 2 ));
1564 double L2 = getDistance(P( 2 ),P( 3 ));
1565 double L3 = getDistance(P( 3 ),P( 4 ));
1566 double L4 = getDistance(P( 4 ),P( 1 ));
1567 double L5 = getDistance(P( 1 ),P( 5 ));
1568 double L6 = getDistance(P( 2 ),P( 5 ));
1569 double L7 = getDistance(P( 3 ),P( 5 ));
1570 double L8 = getDistance(P( 4 ),P( 5 ));
1572 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1573 aVal = Max(aVal,Max(L7,L8));
1576 else if (len == 6){ // pentaidres
1577 double L1 = getDistance(P( 1 ),P( 2 ));
1578 double L2 = getDistance(P( 2 ),P( 3 ));
1579 double L3 = getDistance(P( 3 ),P( 1 ));
1580 double L4 = getDistance(P( 4 ),P( 5 ));
1581 double L5 = getDistance(P( 5 ),P( 6 ));
1582 double L6 = getDistance(P( 6 ),P( 4 ));
1583 double L7 = getDistance(P( 1 ),P( 4 ));
1584 double L8 = getDistance(P( 2 ),P( 5 ));
1585 double L9 = getDistance(P( 3 ),P( 6 ));
1587 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1588 aVal = Max(aVal,Max(Max(L7,L8),L9));
1591 else if (len == 8){ // hexaider
1592 double L1 = getDistance(P( 1 ),P( 2 ));
1593 double L2 = getDistance(P( 2 ),P( 3 ));
1594 double L3 = getDistance(P( 3 ),P( 4 ));
1595 double L4 = getDistance(P( 4 ),P( 1 ));
1596 double L5 = getDistance(P( 5 ),P( 6 ));
1597 double L6 = getDistance(P( 6 ),P( 7 ));
1598 double L7 = getDistance(P( 7 ),P( 8 ));
1599 double L8 = getDistance(P( 8 ),P( 5 ));
1600 double L9 = getDistance(P( 1 ),P( 5 ));
1601 double L10= getDistance(P( 2 ),P( 6 ));
1602 double L11= getDistance(P( 3 ),P( 7 ));
1603 double L12= getDistance(P( 4 ),P( 8 ));
1605 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1606 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1607 aVal = Max(aVal,Max(L11,L12));
1612 if (len == 10){ // quadratic tetraidrs
1613 double L1 = getDistance(P( 1 ),P( 5 )) + getDistance(P( 5 ),P( 2 ));
1614 double L2 = getDistance(P( 2 ),P( 6 )) + getDistance(P( 6 ),P( 3 ));
1615 double L3 = getDistance(P( 3 ),P( 7 )) + getDistance(P( 7 ),P( 1 ));
1616 double L4 = getDistance(P( 1 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1617 double L5 = getDistance(P( 2 ),P( 9 )) + getDistance(P( 9 ),P( 4 ));
1618 double L6 = getDistance(P( 3 ),P( 10 )) + getDistance(P( 10 ),P( 4 ));
1619 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1622 else if (len == 13){ // quadratic piramids
1623 double L1 = getDistance(P( 1 ),P( 6 )) + getDistance(P( 6 ),P( 2 ));
1624 double L2 = getDistance(P( 2 ),P( 7 )) + getDistance(P( 7 ),P( 3 ));
1625 double L3 = getDistance(P( 3 ),P( 8 )) + getDistance(P( 8 ),P( 4 ));
1626 double L4 = getDistance(P( 4 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1627 double L5 = getDistance(P( 1 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1628 double L6 = getDistance(P( 2 ),P( 11 )) + getDistance(P( 11 ),P( 5 ));
1629 double L7 = getDistance(P( 3 ),P( 12 )) + getDistance(P( 12 ),P( 5 ));
1630 double L8 = getDistance(P( 4 ),P( 13 )) + getDistance(P( 13 ),P( 5 ));
1631 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1632 aVal = Max(aVal,Max(L7,L8));
1635 else if (len == 15){ // quadratic pentaidres
1636 double L1 = getDistance(P( 1 ),P( 7 )) + getDistance(P( 7 ),P( 2 ));
1637 double L2 = getDistance(P( 2 ),P( 8 )) + getDistance(P( 8 ),P( 3 ));
1638 double L3 = getDistance(P( 3 ),P( 9 )) + getDistance(P( 9 ),P( 1 ));
1639 double L4 = getDistance(P( 4 ),P( 10 )) + getDistance(P( 10 ),P( 5 ));
1640 double L5 = getDistance(P( 5 ),P( 11 )) + getDistance(P( 11 ),P( 6 ));
1641 double L6 = getDistance(P( 6 ),P( 12 )) + getDistance(P( 12 ),P( 4 ));
1642 double L7 = getDistance(P( 1 ),P( 13 )) + getDistance(P( 13 ),P( 4 ));
1643 double L8 = getDistance(P( 2 ),P( 14 )) + getDistance(P( 14 ),P( 5 ));
1644 double L9 = getDistance(P( 3 ),P( 15 )) + getDistance(P( 15 ),P( 6 ));
1645 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1646 aVal = Max(aVal,Max(Max(L7,L8),L9));
1649 else if (len == 20){ // quadratic hexaider
1650 double L1 = getDistance(P( 1 ),P( 9 )) + getDistance(P( 9 ),P( 2 ));
1651 double L2 = getDistance(P( 2 ),P( 10 )) + getDistance(P( 10 ),P( 3 ));
1652 double L3 = getDistance(P( 3 ),P( 11 )) + getDistance(P( 11 ),P( 4 ));
1653 double L4 = getDistance(P( 4 ),P( 12 )) + getDistance(P( 12 ),P( 1 ));
1654 double L5 = getDistance(P( 5 ),P( 13 )) + getDistance(P( 13 ),P( 6 ));
1655 double L6 = getDistance(P( 6 ),P( 14 )) + getDistance(P( 14 ),P( 7 ));
1656 double L7 = getDistance(P( 7 ),P( 15 )) + getDistance(P( 15 ),P( 8 ));
1657 double L8 = getDistance(P( 8 ),P( 16 )) + getDistance(P( 16 ),P( 5 ));
1658 double L9 = getDistance(P( 1 ),P( 17 )) + getDistance(P( 17 ),P( 5 ));
1659 double L10= getDistance(P( 2 ),P( 18 )) + getDistance(P( 18 ),P( 6 ));
1660 double L11= getDistance(P( 3 ),P( 19 )) + getDistance(P( 19 ),P( 7 ));
1661 double L12= getDistance(P( 4 ),P( 20 )) + getDistance(P( 20 ),P( 8 ));
1662 aVal = Max(Max(Max(L1,L2),Max(L3,L4)),Max(L5,L6));
1663 aVal = Max(aVal,Max(Max(L7,L8),Max(L9,L10)));
1664 aVal = Max(aVal,Max(L11,L12));
1676 if ( myPrecision >= 0 )
1678 double prec = pow( 10., (double)( myPrecision ) );
1679 aVal = floor( aVal * prec + 0.5 ) / prec;
1688 double Length2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1690 // meaningless as it is not quality control functor
1694 SMDSAbs_ElementType Length2D::GetType() const
1696 return SMDSAbs_Face;
1699 Length2D::Value::Value(double theLength,long thePntId1, long thePntId2):
1702 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1703 if(thePntId1 > thePntId2){
1704 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1708 bool Length2D::Value::operator<(const Length2D::Value& x) const{
1709 if(myPntId[0] < x.myPntId[0]) return true;
1710 if(myPntId[0] == x.myPntId[0])
1711 if(myPntId[1] < x.myPntId[1]) return true;
1715 void Length2D::GetValues(TValues& theValues){
1717 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1718 for(; anIter->more(); ){
1719 const SMDS_MeshFace* anElem = anIter->next();
1721 if(anElem->IsQuadratic()) {
1722 const SMDS_VtkFace* F =
1723 dynamic_cast<const SMDS_VtkFace*>(anElem);
1724 // use special nodes iterator
1725 SMDS_ElemIteratorPtr anIter = F->interlacedNodesElemIterator();
1730 const SMDS_MeshElement* aNode;
1732 aNode = anIter->next();
1733 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1734 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1735 aNodeId[0] = aNodeId[1] = aNode->GetID();
1738 for(; anIter->more(); ){
1739 const SMDS_MeshNode* N1 = static_cast<const SMDS_MeshNode*> (anIter->next());
1740 P[2] = gp_Pnt(N1->X(),N1->Y(),N1->Z());
1741 aNodeId[2] = N1->GetID();
1742 aLength = P[1].Distance(P[2]);
1743 if(!anIter->more()) break;
1744 const SMDS_MeshNode* N2 = static_cast<const SMDS_MeshNode*> (anIter->next());
1745 P[3] = gp_Pnt(N2->X(),N2->Y(),N2->Z());
1746 aNodeId[3] = N2->GetID();
1747 aLength += P[2].Distance(P[3]);
1748 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1749 Value aValue2(aLength,aNodeId[2],aNodeId[3]);
1751 aNodeId[1] = aNodeId[3];
1752 theValues.insert(aValue1);
1753 theValues.insert(aValue2);
1755 aLength += P[2].Distance(P[0]);
1756 Value aValue1(aLength,aNodeId[1],aNodeId[2]);
1757 Value aValue2(aLength,aNodeId[2],aNodeId[0]);
1758 theValues.insert(aValue1);
1759 theValues.insert(aValue2);
1762 SMDS_ElemIteratorPtr aNodesIter = anElem->nodesIterator();
1767 const SMDS_MeshElement* aNode;
1768 if(aNodesIter->more()){
1769 aNode = aNodesIter->next();
1770 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1771 P[0] = P[1] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1772 aNodeId[0] = aNodeId[1] = aNode->GetID();
1775 for(; aNodesIter->more(); ){
1776 aNode = aNodesIter->next();
1777 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode;
1778 long anId = aNode->GetID();
1780 P[2] = gp_Pnt(aNodes->X(),aNodes->Y(),aNodes->Z());
1782 aLength = P[1].Distance(P[2]);
1784 Value aValue(aLength,aNodeId[1],anId);
1787 theValues.insert(aValue);
1790 aLength = P[0].Distance(P[1]);
1792 Value aValue(aLength,aNodeId[0],aNodeId[1]);
1793 theValues.insert(aValue);
1799 Class : MultiConnection
1800 Description : Functor for calculating number of faces conneted to the edge
1802 double MultiConnection::GetValue( const TSequenceOfXYZ& P )
1806 double MultiConnection::GetValue( long theId )
1808 return getNbMultiConnection( myMesh, theId );
1811 double MultiConnection::GetBadRate( double Value, int /*nbNodes*/ ) const
1813 // meaningless as it is not quality control functor
1817 SMDSAbs_ElementType MultiConnection::GetType() const
1819 return SMDSAbs_Edge;
1823 Class : MultiConnection2D
1824 Description : Functor for calculating number of faces conneted to the edge
1826 double MultiConnection2D::GetValue( const TSequenceOfXYZ& P )
1831 double MultiConnection2D::GetValue( long theElementId )
1835 const SMDS_MeshElement* aFaceElem = myMesh->FindElement(theElementId);
1836 SMDSAbs_ElementType aType = aFaceElem->GetType();
1841 int i = 0, len = aFaceElem->NbNodes();
1842 SMDS_ElemIteratorPtr anIter = aFaceElem->nodesIterator();
1845 const SMDS_MeshNode *aNode, *aNode0;
1846 TColStd_MapOfInteger aMap, aMapPrev;
1848 for (i = 0; i <= len; i++) {
1853 if (anIter->more()) {
1854 aNode = (SMDS_MeshNode*)anIter->next();
1862 if (i == 0) aNode0 = aNode;
1864 SMDS_ElemIteratorPtr anElemIter = aNode->GetInverseElementIterator();
1865 while (anElemIter->more()) {
1866 const SMDS_MeshElement* anElem = anElemIter->next();
1867 if (anElem != 0 && anElem->GetType() == SMDSAbs_Face) {
1868 int anId = anElem->GetID();
1871 if (aMapPrev.Contains(anId)) {
1876 aResult = Max(aResult, aNb);
1887 double MultiConnection2D::GetBadRate( double Value, int /*nbNodes*/ ) const
1889 // meaningless as it is not quality control functor
1893 SMDSAbs_ElementType MultiConnection2D::GetType() const
1895 return SMDSAbs_Face;
1898 MultiConnection2D::Value::Value(long thePntId1, long thePntId2)
1900 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
1901 if(thePntId1 > thePntId2){
1902 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
1906 bool MultiConnection2D::Value::operator<(const MultiConnection2D::Value& x) const{
1907 if(myPntId[0] < x.myPntId[0]) return true;
1908 if(myPntId[0] == x.myPntId[0])
1909 if(myPntId[1] < x.myPntId[1]) return true;
1913 void MultiConnection2D::GetValues(MValues& theValues){
1914 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
1915 for(; anIter->more(); ){
1916 const SMDS_MeshFace* anElem = anIter->next();
1917 SMDS_ElemIteratorPtr aNodesIter;
1918 if ( anElem->IsQuadratic() )
1919 aNodesIter = dynamic_cast<const SMDS_VtkFace*>
1920 (anElem)->interlacedNodesElemIterator();
1922 aNodesIter = anElem->nodesIterator();
1925 //int aNbConnects=0;
1926 const SMDS_MeshNode* aNode0;
1927 const SMDS_MeshNode* aNode1;
1928 const SMDS_MeshNode* aNode2;
1929 if(aNodesIter->more()){
1930 aNode0 = (SMDS_MeshNode*) aNodesIter->next();
1932 const SMDS_MeshNode* aNodes = (SMDS_MeshNode*) aNode1;
1933 aNodeId[0] = aNodeId[1] = aNodes->GetID();
1935 for(; aNodesIter->more(); ) {
1936 aNode2 = (SMDS_MeshNode*) aNodesIter->next();
1937 long anId = aNode2->GetID();
1940 Value aValue(aNodeId[1],aNodeId[2]);
1941 MValues::iterator aItr = theValues.find(aValue);
1942 if (aItr != theValues.end()){
1947 theValues[aValue] = 1;
1950 //cout << "NodeIds: "<<aNodeId[1]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1951 aNodeId[1] = aNodeId[2];
1954 Value aValue(aNodeId[0],aNodeId[2]);
1955 MValues::iterator aItr = theValues.find(aValue);
1956 if (aItr != theValues.end()) {
1961 theValues[aValue] = 1;
1964 //cout << "NodeIds: "<<aNodeId[0]<<","<<aNodeId[2]<<" nbconn="<<aNbConnects<<endl;
1974 Class : BadOrientedVolume
1975 Description : Predicate bad oriented volumes
1978 BadOrientedVolume::BadOrientedVolume()
1983 void BadOrientedVolume::SetMesh( const SMDS_Mesh* theMesh )
1988 bool BadOrientedVolume::IsSatisfy( long theId )
1993 SMDS_VolumeTool vTool( myMesh->FindElement( theId ));
1994 return !vTool.IsForward();
1997 SMDSAbs_ElementType BadOrientedVolume::GetType() const
1999 return SMDSAbs_Volume;
2003 Class : BareBorderVolume
2006 bool BareBorderVolume::IsSatisfy(long theElementId )
2008 SMDS_VolumeTool myTool;
2009 if ( myTool.Set( myMesh->FindElement(theElementId)))
2011 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2012 if ( myTool.IsFreeFace( iF ))
2014 const SMDS_MeshNode** n = myTool.GetFaceNodes(iF);
2015 vector< const SMDS_MeshNode*> nodes( n, n+myTool.NbFaceNodes(iF));
2016 if ( !myMesh->FindElement( nodes, SMDSAbs_Face, /*Nomedium=*/false))
2024 Class : BareBorderFace
2027 bool BareBorderFace::IsSatisfy(long theElementId )
2030 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2032 if ( face->GetType() == SMDSAbs_Face )
2034 int nbN = face->NbCornerNodes();
2035 for ( int i = 0; i < nbN && !ok; ++i )
2037 // check if a link is shared by another face
2038 const SMDS_MeshNode* n1 = face->GetNode( i );
2039 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2040 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2041 bool isShared = false;
2042 while ( !isShared && fIt->more() )
2044 const SMDS_MeshElement* f = fIt->next();
2045 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2049 const int iQuad = face->IsQuadratic();
2050 myLinkNodes.resize( 2 + iQuad);
2051 myLinkNodes[0] = n1;
2052 myLinkNodes[1] = n2;
2054 myLinkNodes[2] = face->GetNode( i+nbN );
2055 ok = !myMesh->FindElement( myLinkNodes, SMDSAbs_Edge, /*noMedium=*/false);
2064 Class : OverConstrainedVolume
2067 bool OverConstrainedVolume::IsSatisfy(long theElementId )
2069 // An element is over-constrained if it has N-1 free borders where
2070 // N is the number of edges/faces for a 2D/3D element.
2071 SMDS_VolumeTool myTool;
2072 if ( myTool.Set( myMesh->FindElement(theElementId)))
2074 int nbSharedFaces = 0;
2075 for ( int iF = 0; iF < myTool.NbFaces(); ++iF )
2076 if ( !myTool.IsFreeFace( iF ) && ++nbSharedFaces > 1 )
2078 return ( nbSharedFaces == 1 );
2084 Class : OverConstrainedFace
2087 bool OverConstrainedFace::IsSatisfy(long theElementId )
2089 // An element is over-constrained if it has N-1 free borders where
2090 // N is the number of edges/faces for a 2D/3D element.
2091 if ( const SMDS_MeshElement* face = myMesh->FindElement(theElementId))
2092 if ( face->GetType() == SMDSAbs_Face )
2094 int nbSharedBorders = 0;
2095 int nbN = face->NbCornerNodes();
2096 for ( int i = 0; i < nbN; ++i )
2098 // check if a link is shared by another face
2099 const SMDS_MeshNode* n1 = face->GetNode( i );
2100 const SMDS_MeshNode* n2 = face->GetNode( (i+1)%nbN );
2101 SMDS_ElemIteratorPtr fIt = n1->GetInverseElementIterator( SMDSAbs_Face );
2102 bool isShared = false;
2103 while ( !isShared && fIt->more() )
2105 const SMDS_MeshElement* f = fIt->next();
2106 isShared = ( f != face && f->GetNodeIndex(n2) != -1 );
2108 if ( isShared && ++nbSharedBorders > 1 )
2111 return ( nbSharedBorders == 1 );
2117 Class : CoincidentNodes
2118 Description : Predicate of Coincident nodes
2121 CoincidentNodes::CoincidentNodes()
2126 bool CoincidentNodes::IsSatisfy( long theElementId )
2128 return myCoincidentIDs.Contains( theElementId );
2131 SMDSAbs_ElementType CoincidentNodes::GetType() const
2133 return SMDSAbs_Node;
2136 void CoincidentNodes::SetMesh( const SMDS_Mesh* theMesh )
2138 myMeshModifTracer.SetMesh( theMesh );
2139 if ( myMeshModifTracer.IsMeshModified() )
2141 TIDSortedNodeSet nodesToCheck;
2142 SMDS_NodeIteratorPtr nIt = theMesh->nodesIterator(/*idInceasingOrder=*/true);
2143 while ( nIt->more() )
2144 nodesToCheck.insert( nodesToCheck.end(), nIt->next() );
2146 list< list< const SMDS_MeshNode*> > nodeGroups;
2147 SMESH_OctreeNode::FindCoincidentNodes ( nodesToCheck, &nodeGroups, myToler );
2149 myCoincidentIDs.Clear();
2150 list< list< const SMDS_MeshNode*> >::iterator groupIt = nodeGroups.begin();
2151 for ( ; groupIt != nodeGroups.end(); ++groupIt )
2153 list< const SMDS_MeshNode*>& coincNodes = *groupIt;
2154 list< const SMDS_MeshNode*>::iterator n = coincNodes.begin();
2155 for ( ; n != coincNodes.end(); ++n )
2156 myCoincidentIDs.Add( (*n)->GetID() );
2162 Class : CoincidentElements
2163 Description : Predicate of Coincident Elements
2164 Note : This class is suitable only for visualization of Coincident Elements
2167 CoincidentElements::CoincidentElements()
2172 void CoincidentElements::SetMesh( const SMDS_Mesh* theMesh )
2177 bool CoincidentElements::IsSatisfy( long theElementId )
2179 if ( !myMesh ) return false;
2181 if ( const SMDS_MeshElement* e = myMesh->FindElement( theElementId ))
2183 if ( e->GetType() != GetType() ) return false;
2184 set< const SMDS_MeshNode* > elemNodes( e->begin_nodes(), e->end_nodes() );
2185 const int nbNodes = e->NbNodes();
2186 SMDS_ElemIteratorPtr invIt = (*elemNodes.begin())->GetInverseElementIterator( GetType() );
2187 while ( invIt->more() )
2189 const SMDS_MeshElement* e2 = invIt->next();
2190 if ( e2 == e || e2->NbNodes() != nbNodes ) continue;
2192 bool sameNodes = true;
2193 for ( size_t i = 0; i < elemNodes.size() && sameNodes; ++i )
2194 sameNodes = ( elemNodes.count( e2->GetNode( i )));
2202 SMDSAbs_ElementType CoincidentElements1D::GetType() const
2204 return SMDSAbs_Edge;
2206 SMDSAbs_ElementType CoincidentElements2D::GetType() const
2208 return SMDSAbs_Face;
2210 SMDSAbs_ElementType CoincidentElements3D::GetType() const
2212 return SMDSAbs_Volume;
2218 Description : Predicate for free borders
2221 FreeBorders::FreeBorders()
2226 void FreeBorders::SetMesh( const SMDS_Mesh* theMesh )
2231 bool FreeBorders::IsSatisfy( long theId )
2233 return getNbMultiConnection( myMesh, theId ) == 1;
2236 SMDSAbs_ElementType FreeBorders::GetType() const
2238 return SMDSAbs_Edge;
2244 Description : Predicate for free Edges
2246 FreeEdges::FreeEdges()
2251 void FreeEdges::SetMesh( const SMDS_Mesh* theMesh )
2256 bool FreeEdges::IsFreeEdge( const SMDS_MeshNode** theNodes, const int theFaceId )
2258 TColStd_MapOfInteger aMap;
2259 for ( int i = 0; i < 2; i++ )
2261 SMDS_ElemIteratorPtr anElemIter = theNodes[ i ]->GetInverseElementIterator(SMDSAbs_Face);
2262 while( anElemIter->more() )
2264 if ( const SMDS_MeshElement* anElem = anElemIter->next())
2266 const int anId = anElem->GetID();
2267 if ( anId != theFaceId && !aMap.Add( anId ))
2275 bool FreeEdges::IsSatisfy( long theId )
2280 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2281 if ( aFace == 0 || aFace->GetType() != SMDSAbs_Face || aFace->NbNodes() < 3 )
2284 SMDS_ElemIteratorPtr anIter;
2285 if ( aFace->IsQuadratic() ) {
2286 anIter = dynamic_cast<const SMDS_VtkFace*>
2287 (aFace)->interlacedNodesElemIterator();
2290 anIter = aFace->nodesIterator();
2295 int i = 0, nbNodes = aFace->NbNodes();
2296 std::vector <const SMDS_MeshNode*> aNodes( nbNodes+1 );
2297 while( anIter->more() )
2299 const SMDS_MeshNode* aNode = (SMDS_MeshNode*)anIter->next();
2302 aNodes[ i++ ] = aNode;
2304 aNodes[ nbNodes ] = aNodes[ 0 ];
2306 for ( i = 0; i < nbNodes; i++ )
2307 if ( IsFreeEdge( &aNodes[ i ], theId ) )
2313 SMDSAbs_ElementType FreeEdges::GetType() const
2315 return SMDSAbs_Face;
2318 FreeEdges::Border::Border(long theElemId, long thePntId1, long thePntId2):
2321 myPntId[0] = thePntId1; myPntId[1] = thePntId2;
2322 if(thePntId1 > thePntId2){
2323 myPntId[1] = thePntId1; myPntId[0] = thePntId2;
2327 bool FreeEdges::Border::operator<(const FreeEdges::Border& x) const{
2328 if(myPntId[0] < x.myPntId[0]) return true;
2329 if(myPntId[0] == x.myPntId[0])
2330 if(myPntId[1] < x.myPntId[1]) return true;
2334 inline void UpdateBorders(const FreeEdges::Border& theBorder,
2335 FreeEdges::TBorders& theRegistry,
2336 FreeEdges::TBorders& theContainer)
2338 if(theRegistry.find(theBorder) == theRegistry.end()){
2339 theRegistry.insert(theBorder);
2340 theContainer.insert(theBorder);
2342 theContainer.erase(theBorder);
2346 void FreeEdges::GetBoreders(TBorders& theBorders)
2349 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
2350 for(; anIter->more(); ){
2351 const SMDS_MeshFace* anElem = anIter->next();
2352 long anElemId = anElem->GetID();
2353 SMDS_ElemIteratorPtr aNodesIter;
2354 if ( anElem->IsQuadratic() )
2355 aNodesIter = static_cast<const SMDS_VtkFace*>(anElem)->
2356 interlacedNodesElemIterator();
2358 aNodesIter = anElem->nodesIterator();
2360 const SMDS_MeshElement* aNode;
2361 if(aNodesIter->more()){
2362 aNode = aNodesIter->next();
2363 aNodeId[0] = aNodeId[1] = aNode->GetID();
2365 for(; aNodesIter->more(); ){
2366 aNode = aNodesIter->next();
2367 long anId = aNode->GetID();
2368 Border aBorder(anElemId,aNodeId[1],anId);
2370 UpdateBorders(aBorder,aRegistry,theBorders);
2372 Border aBorder(anElemId,aNodeId[0],aNodeId[1]);
2373 UpdateBorders(aBorder,aRegistry,theBorders);
2380 Description : Predicate for free nodes
2383 FreeNodes::FreeNodes()
2388 void FreeNodes::SetMesh( const SMDS_Mesh* theMesh )
2393 bool FreeNodes::IsSatisfy( long theNodeId )
2395 const SMDS_MeshNode* aNode = myMesh->FindNode( theNodeId );
2399 return (aNode->NbInverseElements() < 1);
2402 SMDSAbs_ElementType FreeNodes::GetType() const
2404 return SMDSAbs_Node;
2410 Description : Predicate for free faces
2413 FreeFaces::FreeFaces()
2418 void FreeFaces::SetMesh( const SMDS_Mesh* theMesh )
2423 bool FreeFaces::IsSatisfy( long theId )
2425 if (!myMesh) return false;
2426 // check that faces nodes refers to less than two common volumes
2427 const SMDS_MeshElement* aFace = myMesh->FindElement( theId );
2428 if ( !aFace || aFace->GetType() != SMDSAbs_Face )
2431 int nbNode = aFace->NbNodes();
2433 // collect volumes check that number of volumss with count equal nbNode not less than 2
2434 typedef map< SMDS_MeshElement*, int > TMapOfVolume; // map of volume counters
2435 typedef map< SMDS_MeshElement*, int >::iterator TItrMapOfVolume; // iterator
2436 TMapOfVolume mapOfVol;
2438 SMDS_ElemIteratorPtr nodeItr = aFace->nodesIterator();
2439 while ( nodeItr->more() ) {
2440 const SMDS_MeshNode* aNode = static_cast<const SMDS_MeshNode*>(nodeItr->next());
2441 if ( !aNode ) continue;
2442 SMDS_ElemIteratorPtr volItr = aNode->GetInverseElementIterator(SMDSAbs_Volume);
2443 while ( volItr->more() ) {
2444 SMDS_MeshElement* aVol = (SMDS_MeshElement*)volItr->next();
2445 TItrMapOfVolume itr = mapOfVol.insert(make_pair(aVol, 0)).first;
2450 TItrMapOfVolume volItr = mapOfVol.begin();
2451 TItrMapOfVolume volEnd = mapOfVol.end();
2452 for ( ; volItr != volEnd; ++volItr )
2453 if ( (*volItr).second >= nbNode )
2455 // face is not free if number of volumes constructed on thier nodes more than one
2459 SMDSAbs_ElementType FreeFaces::GetType() const
2461 return SMDSAbs_Face;
2465 Class : LinearOrQuadratic
2466 Description : Predicate to verify whether a mesh element is linear
2469 LinearOrQuadratic::LinearOrQuadratic()
2474 void LinearOrQuadratic::SetMesh( const SMDS_Mesh* theMesh )
2479 bool LinearOrQuadratic::IsSatisfy( long theId )
2481 if (!myMesh) return false;
2482 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2483 if ( !anElem || (myType != SMDSAbs_All && anElem->GetType() != myType) )
2485 return (!anElem->IsQuadratic());
2488 void LinearOrQuadratic::SetType( SMDSAbs_ElementType theType )
2493 SMDSAbs_ElementType LinearOrQuadratic::GetType() const
2500 Description : Functor for check color of group to whic mesh element belongs to
2503 GroupColor::GroupColor()
2507 bool GroupColor::IsSatisfy( long theId )
2509 return (myIDs.find( theId ) != myIDs.end());
2512 void GroupColor::SetType( SMDSAbs_ElementType theType )
2517 SMDSAbs_ElementType GroupColor::GetType() const
2522 static bool isEqual( const Quantity_Color& theColor1,
2523 const Quantity_Color& theColor2 )
2525 // tolerance to compare colors
2526 const double tol = 5*1e-3;
2527 return ( fabs( theColor1.Red() - theColor2.Red() ) < tol &&
2528 fabs( theColor1.Green() - theColor2.Green() ) < tol &&
2529 fabs( theColor1.Blue() - theColor2.Blue() ) < tol );
2533 void GroupColor::SetMesh( const SMDS_Mesh* theMesh )
2537 const SMESHDS_Mesh* aMesh = dynamic_cast<const SMESHDS_Mesh*>(theMesh);
2541 int nbGrp = aMesh->GetNbGroups();
2545 // iterates on groups and find necessary elements ids
2546 const std::set<SMESHDS_GroupBase*>& aGroups = aMesh->GetGroups();
2547 set<SMESHDS_GroupBase*>::const_iterator GrIt = aGroups.begin();
2548 for (; GrIt != aGroups.end(); GrIt++) {
2549 SMESHDS_GroupBase* aGrp = (*GrIt);
2552 // check type and color of group
2553 if ( !isEqual( myColor, aGrp->GetColor() ) )
2555 if ( myType != SMDSAbs_All && myType != (SMDSAbs_ElementType)aGrp->GetType() )
2558 SMDSAbs_ElementType aGrpElType = (SMDSAbs_ElementType)aGrp->GetType();
2559 if ( myType == aGrpElType || (myType == SMDSAbs_All && aGrpElType != SMDSAbs_Node) ) {
2560 // add elements IDS into control
2561 int aSize = aGrp->Extent();
2562 for (int i = 0; i < aSize; i++)
2563 myIDs.insert( aGrp->GetID(i+1) );
2568 void GroupColor::SetColorStr( const TCollection_AsciiString& theStr )
2570 TCollection_AsciiString aStr = theStr;
2571 aStr.RemoveAll( ' ' );
2572 aStr.RemoveAll( '\t' );
2573 for ( int aPos = aStr.Search( ";;" ); aPos != -1; aPos = aStr.Search( ";;" ) )
2574 aStr.Remove( aPos, 2 );
2575 Standard_Real clr[3];
2576 clr[0] = clr[1] = clr[2] = 0.;
2577 for ( int i = 0; i < 3; i++ ) {
2578 TCollection_AsciiString tmpStr = aStr.Token( ";", i+1 );
2579 if ( !tmpStr.IsEmpty() && tmpStr.IsRealValue() )
2580 clr[i] = tmpStr.RealValue();
2582 myColor = Quantity_Color( clr[0], clr[1], clr[2], Quantity_TOC_RGB );
2585 //=======================================================================
2586 // name : GetRangeStr
2587 // Purpose : Get range as a string.
2588 // Example: "1,2,3,50-60,63,67,70-"
2589 //=======================================================================
2590 void GroupColor::GetColorStr( TCollection_AsciiString& theResStr ) const
2593 theResStr += TCollection_AsciiString( myColor.Red() );
2594 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Green() );
2595 theResStr += TCollection_AsciiString( ";" ) + TCollection_AsciiString( myColor.Blue() );
2599 Class : ElemGeomType
2600 Description : Predicate to check element geometry type
2603 ElemGeomType::ElemGeomType()
2606 myType = SMDSAbs_All;
2607 myGeomType = SMDSGeom_TRIANGLE;
2610 void ElemGeomType::SetMesh( const SMDS_Mesh* theMesh )
2615 bool ElemGeomType::IsSatisfy( long theId )
2617 if (!myMesh) return false;
2618 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2621 const SMDSAbs_ElementType anElemType = anElem->GetType();
2622 if ( myType != SMDSAbs_All && anElemType != myType )
2624 const int aNbNode = anElem->NbNodes();
2626 switch( anElemType )
2629 isOk = (myGeomType == SMDSGeom_POINT);
2633 isOk = (myGeomType == SMDSGeom_EDGE);
2637 if ( myGeomType == SMDSGeom_TRIANGLE )
2638 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 6 : aNbNode == 3));
2639 else if ( myGeomType == SMDSGeom_QUADRANGLE )
2640 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? ( aNbNode == 8 || aNbNode == 9 ) : aNbNode == 4));
2641 else if ( myGeomType == SMDSGeom_POLYGON )
2642 isOk = anElem->IsPoly();
2645 case SMDSAbs_Volume:
2646 if ( myGeomType == SMDSGeom_TETRA )
2647 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 10 : aNbNode == 4));
2648 else if ( myGeomType == SMDSGeom_PYRAMID )
2649 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 13 : aNbNode == 5));
2650 else if ( myGeomType == SMDSGeom_PENTA )
2651 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? aNbNode == 15 : aNbNode == 6));
2652 else if ( myGeomType == SMDSGeom_HEXA )
2653 isOk = (!anElem->IsPoly() && (anElem->IsQuadratic() ? ( aNbNode == 20 || aNbNode == 27 ): aNbNode == 8));
2654 else if ( myGeomType == SMDSGeom_HEXAGONAL_PRISM )
2655 isOk = (anElem->GetEntityType() == SMDSEntity_Hexagonal_Prism );
2656 else if ( myGeomType == SMDSGeom_POLYHEDRA )
2657 isOk = anElem->IsPoly();
2664 void ElemGeomType::SetType( SMDSAbs_ElementType theType )
2669 SMDSAbs_ElementType ElemGeomType::GetType() const
2674 void ElemGeomType::SetGeomType( SMDSAbs_GeometryType theType )
2676 myGeomType = theType;
2679 SMDSAbs_GeometryType ElemGeomType::GetGeomType() const
2684 //================================================================================
2686 * \brief Class CoplanarFaces
2688 //================================================================================
2690 CoplanarFaces::CoplanarFaces()
2691 : myFaceID(0), myToler(0)
2694 void CoplanarFaces::SetMesh( const SMDS_Mesh* theMesh )
2696 myMeshModifTracer.SetMesh( theMesh );
2697 if ( myMeshModifTracer.IsMeshModified() )
2699 // Build a set of coplanar face ids
2701 myCoplanarIDs.clear();
2703 if ( !myMeshModifTracer.GetMesh() || !myFaceID || !myToler )
2706 const SMDS_MeshElement* face = myMeshModifTracer.GetMesh()->FindElement( myFaceID );
2707 if ( !face || face->GetType() != SMDSAbs_Face )
2711 gp_Vec myNorm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
2715 const double radianTol = myToler * M_PI / 180.;
2716 typedef SMDS_StdIterator< const SMDS_MeshElement*, SMDS_ElemIteratorPtr > TFaceIt;
2717 std::set<const SMDS_MeshElement*> checkedFaces, checkedNodes;
2718 std::list<const SMDS_MeshElement*> faceQueue( 1, face );
2719 while ( !faceQueue.empty() )
2721 face = faceQueue.front();
2722 if ( checkedFaces.insert( face ).second )
2724 gp_Vec norm = getNormale( static_cast<const SMDS_MeshFace*>(face), &normOK );
2725 if (!normOK || myNorm.Angle( norm ) <= radianTol)
2727 myCoplanarIDs.insert( face->GetID() );
2728 std::set<const SMDS_MeshElement*> neighborFaces;
2729 for ( int i = 0; i < face->NbCornerNodes(); ++i )
2731 const SMDS_MeshNode* n = face->GetNode( i );
2732 if ( checkedNodes.insert( n ).second )
2733 neighborFaces.insert( TFaceIt( n->GetInverseElementIterator(SMDSAbs_Face)),
2736 faceQueue.insert( faceQueue.end(), neighborFaces.begin(), neighborFaces.end() );
2739 faceQueue.pop_front();
2743 bool CoplanarFaces::IsSatisfy( long theElementId )
2745 return myCoplanarIDs.count( theElementId );
2750 *Description : Predicate for Range of Ids.
2751 * Range may be specified with two ways.
2752 * 1. Using AddToRange method
2753 * 2. With SetRangeStr method. Parameter of this method is a string
2754 * like as "1,2,3,50-60,63,67,70-"
2757 //=======================================================================
2758 // name : RangeOfIds
2759 // Purpose : Constructor
2760 //=======================================================================
2761 RangeOfIds::RangeOfIds()
2764 myType = SMDSAbs_All;
2767 //=======================================================================
2769 // Purpose : Set mesh
2770 //=======================================================================
2771 void RangeOfIds::SetMesh( const SMDS_Mesh* theMesh )
2776 //=======================================================================
2777 // name : AddToRange
2778 // Purpose : Add ID to the range
2779 //=======================================================================
2780 bool RangeOfIds::AddToRange( long theEntityId )
2782 myIds.Add( theEntityId );
2786 //=======================================================================
2787 // name : GetRangeStr
2788 // Purpose : Get range as a string.
2789 // Example: "1,2,3,50-60,63,67,70-"
2790 //=======================================================================
2791 void RangeOfIds::GetRangeStr( TCollection_AsciiString& theResStr )
2795 TColStd_SequenceOfInteger anIntSeq;
2796 TColStd_SequenceOfAsciiString aStrSeq;
2798 TColStd_MapIteratorOfMapOfInteger anIter( myIds );
2799 for ( ; anIter.More(); anIter.Next() )
2801 int anId = anIter.Key();
2802 TCollection_AsciiString aStr( anId );
2803 anIntSeq.Append( anId );
2804 aStrSeq.Append( aStr );
2807 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2809 int aMinId = myMin( i );
2810 int aMaxId = myMax( i );
2812 TCollection_AsciiString aStr;
2813 if ( aMinId != IntegerFirst() )
2818 if ( aMaxId != IntegerLast() )
2821 // find position of the string in result sequence and insert string in it
2822 if ( anIntSeq.Length() == 0 )
2824 anIntSeq.Append( aMinId );
2825 aStrSeq.Append( aStr );
2829 if ( aMinId < anIntSeq.First() )
2831 anIntSeq.Prepend( aMinId );
2832 aStrSeq.Prepend( aStr );
2834 else if ( aMinId > anIntSeq.Last() )
2836 anIntSeq.Append( aMinId );
2837 aStrSeq.Append( aStr );
2840 for ( int j = 1, k = anIntSeq.Length(); j <= k; j++ )
2841 if ( aMinId < anIntSeq( j ) )
2843 anIntSeq.InsertBefore( j, aMinId );
2844 aStrSeq.InsertBefore( j, aStr );
2850 if ( aStrSeq.Length() == 0 )
2853 theResStr = aStrSeq( 1 );
2854 for ( int j = 2, k = aStrSeq.Length(); j <= k; j++ )
2857 theResStr += aStrSeq( j );
2861 //=======================================================================
2862 // name : SetRangeStr
2863 // Purpose : Define range with string
2864 // Example of entry string: "1,2,3,50-60,63,67,70-"
2865 //=======================================================================
2866 bool RangeOfIds::SetRangeStr( const TCollection_AsciiString& theStr )
2872 TCollection_AsciiString aStr = theStr;
2873 aStr.RemoveAll( ' ' );
2874 aStr.RemoveAll( '\t' );
2876 for ( int aPos = aStr.Search( ",," ); aPos != -1; aPos = aStr.Search( ",," ) )
2877 aStr.Remove( aPos, 2 );
2879 TCollection_AsciiString tmpStr = aStr.Token( ",", 1 );
2881 while ( tmpStr != "" )
2883 tmpStr = aStr.Token( ",", i++ );
2884 int aPos = tmpStr.Search( '-' );
2888 if ( tmpStr.IsIntegerValue() )
2889 myIds.Add( tmpStr.IntegerValue() );
2895 TCollection_AsciiString aMaxStr = tmpStr.Split( aPos );
2896 TCollection_AsciiString aMinStr = tmpStr;
2898 while ( aMinStr.Search( "-" ) != -1 ) aMinStr.RemoveAll( '-' );
2899 while ( aMaxStr.Search( "-" ) != -1 ) aMaxStr.RemoveAll( '-' );
2901 if ( (!aMinStr.IsEmpty() && !aMinStr.IsIntegerValue()) ||
2902 (!aMaxStr.IsEmpty() && !aMaxStr.IsIntegerValue()) )
2905 myMin.Append( aMinStr.IsEmpty() ? IntegerFirst() : aMinStr.IntegerValue() );
2906 myMax.Append( aMaxStr.IsEmpty() ? IntegerLast() : aMaxStr.IntegerValue() );
2913 //=======================================================================
2915 // Purpose : Get type of supported entities
2916 //=======================================================================
2917 SMDSAbs_ElementType RangeOfIds::GetType() const
2922 //=======================================================================
2924 // Purpose : Set type of supported entities
2925 //=======================================================================
2926 void RangeOfIds::SetType( SMDSAbs_ElementType theType )
2931 //=======================================================================
2933 // Purpose : Verify whether entity satisfies to this rpedicate
2934 //=======================================================================
2935 bool RangeOfIds::IsSatisfy( long theId )
2940 if ( myType == SMDSAbs_Node )
2942 if ( myMesh->FindNode( theId ) == 0 )
2947 const SMDS_MeshElement* anElem = myMesh->FindElement( theId );
2948 if ( anElem == 0 || (myType != anElem->GetType() && myType != SMDSAbs_All ))
2952 if ( myIds.Contains( theId ) )
2955 for ( int i = 1, n = myMin.Length(); i <= n; i++ )
2956 if ( theId >= myMin( i ) && theId <= myMax( i ) )
2964 Description : Base class for comparators
2966 Comparator::Comparator():
2970 Comparator::~Comparator()
2973 void Comparator::SetMesh( const SMDS_Mesh* theMesh )
2976 myFunctor->SetMesh( theMesh );
2979 void Comparator::SetMargin( double theValue )
2981 myMargin = theValue;
2984 void Comparator::SetNumFunctor( NumericalFunctorPtr theFunct )
2986 myFunctor = theFunct;
2989 SMDSAbs_ElementType Comparator::GetType() const
2991 return myFunctor ? myFunctor->GetType() : SMDSAbs_All;
2994 double Comparator::GetMargin()
3002 Description : Comparator "<"
3004 bool LessThan::IsSatisfy( long theId )
3006 return myFunctor && myFunctor->GetValue( theId ) < myMargin;
3012 Description : Comparator ">"
3014 bool MoreThan::IsSatisfy( long theId )
3016 return myFunctor && myFunctor->GetValue( theId ) > myMargin;
3022 Description : Comparator "="
3025 myToler(Precision::Confusion())
3028 bool EqualTo::IsSatisfy( long theId )
3030 return myFunctor && fabs( myFunctor->GetValue( theId ) - myMargin ) < myToler;
3033 void EqualTo::SetTolerance( double theToler )
3038 double EqualTo::GetTolerance()
3045 Description : Logical NOT predicate
3047 LogicalNOT::LogicalNOT()
3050 LogicalNOT::~LogicalNOT()
3053 bool LogicalNOT::IsSatisfy( long theId )
3055 return myPredicate && !myPredicate->IsSatisfy( theId );
3058 void LogicalNOT::SetMesh( const SMDS_Mesh* theMesh )
3061 myPredicate->SetMesh( theMesh );
3064 void LogicalNOT::SetPredicate( PredicatePtr thePred )
3066 myPredicate = thePred;
3069 SMDSAbs_ElementType LogicalNOT::GetType() const
3071 return myPredicate ? myPredicate->GetType() : SMDSAbs_All;
3076 Class : LogicalBinary
3077 Description : Base class for binary logical predicate
3079 LogicalBinary::LogicalBinary()
3082 LogicalBinary::~LogicalBinary()
3085 void LogicalBinary::SetMesh( const SMDS_Mesh* theMesh )
3088 myPredicate1->SetMesh( theMesh );
3091 myPredicate2->SetMesh( theMesh );
3094 void LogicalBinary::SetPredicate1( PredicatePtr thePredicate )
3096 myPredicate1 = thePredicate;
3099 void LogicalBinary::SetPredicate2( PredicatePtr thePredicate )
3101 myPredicate2 = thePredicate;
3104 SMDSAbs_ElementType LogicalBinary::GetType() const
3106 if ( !myPredicate1 || !myPredicate2 )
3109 SMDSAbs_ElementType aType1 = myPredicate1->GetType();
3110 SMDSAbs_ElementType aType2 = myPredicate2->GetType();
3112 return aType1 == aType2 ? aType1 : SMDSAbs_All;
3118 Description : Logical AND
3120 bool LogicalAND::IsSatisfy( long theId )
3125 myPredicate1->IsSatisfy( theId ) &&
3126 myPredicate2->IsSatisfy( theId );
3132 Description : Logical OR
3134 bool LogicalOR::IsSatisfy( long theId )
3139 (myPredicate1->IsSatisfy( theId ) ||
3140 myPredicate2->IsSatisfy( theId ));
3154 void Filter::SetPredicate( PredicatePtr thePredicate )
3156 myPredicate = thePredicate;
3159 template<class TElement, class TIterator, class TPredicate>
3160 inline void FillSequence(const TIterator& theIterator,
3161 TPredicate& thePredicate,
3162 Filter::TIdSequence& theSequence)
3164 if ( theIterator ) {
3165 while( theIterator->more() ) {
3166 TElement anElem = theIterator->next();
3167 long anId = anElem->GetID();
3168 if ( thePredicate->IsSatisfy( anId ) )
3169 theSequence.push_back( anId );
3176 GetElementsId( const SMDS_Mesh* theMesh,
3177 PredicatePtr thePredicate,
3178 TIdSequence& theSequence )
3180 theSequence.clear();
3182 if ( !theMesh || !thePredicate )
3185 thePredicate->SetMesh( theMesh );
3187 SMDSAbs_ElementType aType = thePredicate->GetType();
3190 FillSequence<const SMDS_MeshNode*>(theMesh->nodesIterator(),thePredicate,theSequence);
3193 FillSequence<const SMDS_MeshElement*>(theMesh->edgesIterator(),thePredicate,theSequence);
3196 FillSequence<const SMDS_MeshElement*>(theMesh->facesIterator(),thePredicate,theSequence);
3198 case SMDSAbs_Volume:
3199 FillSequence<const SMDS_MeshElement*>(theMesh->volumesIterator(),thePredicate,theSequence);
3202 FillSequence<const SMDS_MeshElement*>(theMesh->edgesIterator(),thePredicate,theSequence);
3203 FillSequence<const SMDS_MeshElement*>(theMesh->facesIterator(),thePredicate,theSequence);
3204 FillSequence<const SMDS_MeshElement*>(theMesh->volumesIterator(),thePredicate,theSequence);
3210 Filter::GetElementsId( const SMDS_Mesh* theMesh,
3211 Filter::TIdSequence& theSequence )
3213 GetElementsId(theMesh,myPredicate,theSequence);
3220 typedef std::set<SMDS_MeshFace*> TMapOfFacePtr;
3226 ManifoldPart::Link::Link( SMDS_MeshNode* theNode1,
3227 SMDS_MeshNode* theNode2 )
3233 ManifoldPart::Link::~Link()
3239 bool ManifoldPart::Link::IsEqual( const ManifoldPart::Link& theLink ) const
3241 if ( myNode1 == theLink.myNode1 &&
3242 myNode2 == theLink.myNode2 )
3244 else if ( myNode1 == theLink.myNode2 &&
3245 myNode2 == theLink.myNode1 )
3251 bool ManifoldPart::Link::operator<( const ManifoldPart::Link& x ) const
3253 if(myNode1 < x.myNode1) return true;
3254 if(myNode1 == x.myNode1)
3255 if(myNode2 < x.myNode2) return true;
3259 bool ManifoldPart::IsEqual( const ManifoldPart::Link& theLink1,
3260 const ManifoldPart::Link& theLink2 )
3262 return theLink1.IsEqual( theLink2 );
3265 ManifoldPart::ManifoldPart()
3268 myAngToler = Precision::Angular();
3269 myIsOnlyManifold = true;
3272 ManifoldPart::~ManifoldPart()
3277 void ManifoldPart::SetMesh( const SMDS_Mesh* theMesh )
3283 SMDSAbs_ElementType ManifoldPart::GetType() const
3284 { return SMDSAbs_Face; }
3286 bool ManifoldPart::IsSatisfy( long theElementId )
3288 return myMapIds.Contains( theElementId );
3291 void ManifoldPart::SetAngleTolerance( const double theAngToler )
3292 { myAngToler = theAngToler; }
3294 double ManifoldPart::GetAngleTolerance() const
3295 { return myAngToler; }
3297 void ManifoldPart::SetIsOnlyManifold( const bool theIsOnly )
3298 { myIsOnlyManifold = theIsOnly; }
3300 void ManifoldPart::SetStartElem( const long theStartId )
3301 { myStartElemId = theStartId; }
3303 bool ManifoldPart::process()
3306 myMapBadGeomIds.Clear();
3308 myAllFacePtr.clear();
3309 myAllFacePtrIntDMap.clear();
3313 // collect all faces into own map
3314 SMDS_FaceIteratorPtr anFaceItr = myMesh->facesIterator();
3315 for (; anFaceItr->more(); )
3317 SMDS_MeshFace* aFacePtr = (SMDS_MeshFace*)anFaceItr->next();
3318 myAllFacePtr.push_back( aFacePtr );
3319 myAllFacePtrIntDMap[aFacePtr] = myAllFacePtr.size()-1;
3322 SMDS_MeshFace* aStartFace = (SMDS_MeshFace*)myMesh->FindElement( myStartElemId );
3326 // the map of non manifold links and bad geometry
3327 TMapOfLink aMapOfNonManifold;
3328 TColStd_MapOfInteger aMapOfTreated;
3330 // begin cycle on faces from start index and run on vector till the end
3331 // and from begin to start index to cover whole vector
3332 const int aStartIndx = myAllFacePtrIntDMap[aStartFace];
3333 bool isStartTreat = false;
3334 for ( int fi = aStartIndx; !isStartTreat || fi != aStartIndx ; fi++ )
3336 if ( fi == aStartIndx )
3337 isStartTreat = true;
3338 // as result next time when fi will be equal to aStartIndx
3340 SMDS_MeshFace* aFacePtr = myAllFacePtr[ fi ];
3341 if ( aMapOfTreated.Contains( aFacePtr->GetID() ) )
3344 aMapOfTreated.Add( aFacePtr->GetID() );
3345 TColStd_MapOfInteger aResFaces;
3346 if ( !findConnected( myAllFacePtrIntDMap, aFacePtr,
3347 aMapOfNonManifold, aResFaces ) )
3349 TColStd_MapIteratorOfMapOfInteger anItr( aResFaces );
3350 for ( ; anItr.More(); anItr.Next() )
3352 int aFaceId = anItr.Key();
3353 aMapOfTreated.Add( aFaceId );
3354 myMapIds.Add( aFaceId );
3357 if ( fi == ( myAllFacePtr.size() - 1 ) )
3359 } // end run on vector of faces
3360 return !myMapIds.IsEmpty();
3363 static void getLinks( const SMDS_MeshFace* theFace,
3364 ManifoldPart::TVectorOfLink& theLinks )
3366 int aNbNode = theFace->NbNodes();
3367 SMDS_ElemIteratorPtr aNodeItr = theFace->nodesIterator();
3369 SMDS_MeshNode* aNode = 0;
3370 for ( ; aNodeItr->more() && i <= aNbNode; )
3373 SMDS_MeshNode* aN1 = (SMDS_MeshNode*)aNodeItr->next();
3377 SMDS_MeshNode* aN2 = ( i >= aNbNode ) ? aNode : (SMDS_MeshNode*)aNodeItr->next();
3379 ManifoldPart::Link aLink( aN1, aN2 );
3380 theLinks.push_back( aLink );
3384 bool ManifoldPart::findConnected
3385 ( const ManifoldPart::TDataMapFacePtrInt& theAllFacePtrInt,
3386 SMDS_MeshFace* theStartFace,
3387 ManifoldPart::TMapOfLink& theNonManifold,
3388 TColStd_MapOfInteger& theResFaces )
3390 theResFaces.Clear();
3391 if ( !theAllFacePtrInt.size() )
3394 if ( getNormale( theStartFace ).SquareModulus() <= gp::Resolution() )
3396 myMapBadGeomIds.Add( theStartFace->GetID() );
3400 ManifoldPart::TMapOfLink aMapOfBoundary, aMapToSkip;
3401 ManifoldPart::TVectorOfLink aSeqOfBoundary;
3402 theResFaces.Add( theStartFace->GetID() );
3403 ManifoldPart::TDataMapOfLinkFacePtr aDMapLinkFace;
3405 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3406 aDMapLinkFace, theNonManifold, theStartFace );
3408 bool isDone = false;
3409 while ( !isDone && aMapOfBoundary.size() != 0 )
3411 bool isToReset = false;
3412 ManifoldPart::TVectorOfLink::iterator pLink = aSeqOfBoundary.begin();
3413 for ( ; !isToReset && pLink != aSeqOfBoundary.end(); ++pLink )
3415 ManifoldPart::Link aLink = *pLink;
3416 if ( aMapToSkip.find( aLink ) != aMapToSkip.end() )
3418 // each link could be treated only once
3419 aMapToSkip.insert( aLink );
3421 ManifoldPart::TVectorOfFacePtr aFaces;
3423 if ( myIsOnlyManifold &&
3424 (theNonManifold.find( aLink ) != theNonManifold.end()) )
3428 getFacesByLink( aLink, aFaces );
3429 // filter the element to keep only indicated elements
3430 ManifoldPart::TVectorOfFacePtr aFiltered;
3431 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3432 for ( ; pFace != aFaces.end(); ++pFace )
3434 SMDS_MeshFace* aFace = *pFace;
3435 if ( myAllFacePtrIntDMap.find( aFace ) != myAllFacePtrIntDMap.end() )
3436 aFiltered.push_back( aFace );
3439 if ( aFaces.size() < 2 ) // no neihgbour faces
3441 else if ( myIsOnlyManifold && aFaces.size() > 2 ) // non manifold case
3443 theNonManifold.insert( aLink );
3448 // compare normal with normals of neighbor element
3449 SMDS_MeshFace* aPrevFace = aDMapLinkFace[ aLink ];
3450 ManifoldPart::TVectorOfFacePtr::iterator pFace = aFaces.begin();
3451 for ( ; pFace != aFaces.end(); ++pFace )
3453 SMDS_MeshFace* aNextFace = *pFace;
3454 if ( aPrevFace == aNextFace )
3456 int anNextFaceID = aNextFace->GetID();
3457 if ( myIsOnlyManifold && theResFaces.Contains( anNextFaceID ) )
3458 // should not be with non manifold restriction. probably bad topology
3460 // check if face was treated and skipped
3461 if ( myMapBadGeomIds.Contains( anNextFaceID ) ||
3462 !isInPlane( aPrevFace, aNextFace ) )
3464 // add new element to connected and extend the boundaries.
3465 theResFaces.Add( anNextFaceID );
3466 expandBoundary( aMapOfBoundary, aSeqOfBoundary,
3467 aDMapLinkFace, theNonManifold, aNextFace );
3471 isDone = !isToReset;
3474 return !theResFaces.IsEmpty();
3477 bool ManifoldPart::isInPlane( const SMDS_MeshFace* theFace1,
3478 const SMDS_MeshFace* theFace2 )
3480 gp_Dir aNorm1 = gp_Dir( getNormale( theFace1 ) );
3481 gp_XYZ aNorm2XYZ = getNormale( theFace2 );
3482 if ( aNorm2XYZ.SquareModulus() <= gp::Resolution() )
3484 myMapBadGeomIds.Add( theFace2->GetID() );
3487 if ( aNorm1.IsParallel( gp_Dir( aNorm2XYZ ), myAngToler ) )
3493 void ManifoldPart::expandBoundary
3494 ( ManifoldPart::TMapOfLink& theMapOfBoundary,
3495 ManifoldPart::TVectorOfLink& theSeqOfBoundary,
3496 ManifoldPart::TDataMapOfLinkFacePtr& theDMapLinkFacePtr,
3497 ManifoldPart::TMapOfLink& theNonManifold,
3498 SMDS_MeshFace* theNextFace ) const
3500 ManifoldPart::TVectorOfLink aLinks;
3501 getLinks( theNextFace, aLinks );
3502 int aNbLink = (int)aLinks.size();
3503 for ( int i = 0; i < aNbLink; i++ )
3505 ManifoldPart::Link aLink = aLinks[ i ];
3506 if ( myIsOnlyManifold && (theNonManifold.find( aLink ) != theNonManifold.end()) )
3508 if ( theMapOfBoundary.find( aLink ) != theMapOfBoundary.end() )
3510 if ( myIsOnlyManifold )
3512 // remove from boundary
3513 theMapOfBoundary.erase( aLink );
3514 ManifoldPart::TVectorOfLink::iterator pLink = theSeqOfBoundary.begin();
3515 for ( ; pLink != theSeqOfBoundary.end(); ++pLink )
3517 ManifoldPart::Link aBoundLink = *pLink;
3518 if ( aBoundLink.IsEqual( aLink ) )
3520 theSeqOfBoundary.erase( pLink );
3528 theMapOfBoundary.insert( aLink );
3529 theSeqOfBoundary.push_back( aLink );
3530 theDMapLinkFacePtr[ aLink ] = theNextFace;
3535 void ManifoldPart::getFacesByLink( const ManifoldPart::Link& theLink,
3536 ManifoldPart::TVectorOfFacePtr& theFaces ) const
3538 std::set<SMDS_MeshCell *> aSetOfFaces;
3539 // take all faces that shared first node
3540 SMDS_ElemIteratorPtr anItr = theLink.myNode1->facesIterator();
3541 for ( ; anItr->more(); )
3543 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3546 aSetOfFaces.insert( aFace );
3548 // take all faces that shared second node
3549 anItr = theLink.myNode2->facesIterator();
3550 // find the common part of two sets
3551 for ( ; anItr->more(); )
3553 SMDS_MeshFace* aFace = (SMDS_MeshFace*)anItr->next();
3554 if ( aSetOfFaces.count( aFace ) )
3555 theFaces.push_back( aFace );
3564 ElementsOnSurface::ElementsOnSurface()
3568 myType = SMDSAbs_All;
3570 myToler = Precision::Confusion();
3571 myUseBoundaries = false;
3574 ElementsOnSurface::~ElementsOnSurface()
3579 void ElementsOnSurface::SetMesh( const SMDS_Mesh* theMesh )
3581 if ( myMesh == theMesh )
3587 bool ElementsOnSurface::IsSatisfy( long theElementId )
3589 return myIds.Contains( theElementId );
3592 SMDSAbs_ElementType ElementsOnSurface::GetType() const
3595 void ElementsOnSurface::SetTolerance( const double theToler )
3597 if ( myToler != theToler )
3602 double ElementsOnSurface::GetTolerance() const
3605 void ElementsOnSurface::SetUseBoundaries( bool theUse )
3607 if ( myUseBoundaries != theUse ) {
3608 myUseBoundaries = theUse;
3609 SetSurface( mySurf, myType );
3613 void ElementsOnSurface::SetSurface( const TopoDS_Shape& theShape,
3614 const SMDSAbs_ElementType theType )
3619 if ( theShape.IsNull() || theShape.ShapeType() != TopAbs_FACE )
3621 mySurf = TopoDS::Face( theShape );
3622 BRepAdaptor_Surface SA( mySurf, myUseBoundaries );
3624 u1 = SA.FirstUParameter(),
3625 u2 = SA.LastUParameter(),
3626 v1 = SA.FirstVParameter(),
3627 v2 = SA.LastVParameter();
3628 Handle(Geom_Surface) surf = BRep_Tool::Surface( mySurf );
3629 myProjector.Init( surf, u1,u2, v1,v2 );
3633 void ElementsOnSurface::process()
3636 if ( mySurf.IsNull() )
3642 if ( myType == SMDSAbs_Face || myType == SMDSAbs_All )
3644 myIds.ReSize( myMesh->NbFaces() );
3645 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
3646 for(; anIter->more(); )
3647 process( anIter->next() );
3650 if ( myType == SMDSAbs_Edge || myType == SMDSAbs_All )
3652 myIds.ReSize( myIds.Extent() + myMesh->NbEdges() );
3653 SMDS_EdgeIteratorPtr anIter = myMesh->edgesIterator();
3654 for(; anIter->more(); )
3655 process( anIter->next() );
3658 if ( myType == SMDSAbs_Node )
3660 myIds.ReSize( myMesh->NbNodes() );
3661 SMDS_NodeIteratorPtr anIter = myMesh->nodesIterator();
3662 for(; anIter->more(); )
3663 process( anIter->next() );
3667 void ElementsOnSurface::process( const SMDS_MeshElement* theElemPtr )
3669 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3670 bool isSatisfy = true;
3671 for ( ; aNodeItr->more(); )
3673 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3674 if ( !isOnSurface( aNode ) )
3681 myIds.Add( theElemPtr->GetID() );
3684 bool ElementsOnSurface::isOnSurface( const SMDS_MeshNode* theNode )
3686 if ( mySurf.IsNull() )
3689 gp_Pnt aPnt( theNode->X(), theNode->Y(), theNode->Z() );
3690 // double aToler2 = myToler * myToler;
3691 // if ( mySurf->IsKind(STANDARD_TYPE(Geom_Plane)))
3693 // gp_Pln aPln = Handle(Geom_Plane)::DownCast(mySurf)->Pln();
3694 // if ( aPln.SquareDistance( aPnt ) > aToler2 )
3697 // else if ( mySurf->IsKind(STANDARD_TYPE(Geom_CylindricalSurface)))
3699 // gp_Cylinder aCyl = Handle(Geom_CylindricalSurface)::DownCast(mySurf)->Cylinder();
3700 // double aRad = aCyl.Radius();
3701 // gp_Ax3 anAxis = aCyl.Position();
3702 // gp_XYZ aLoc = aCyl.Location().XYZ();
3703 // double aXDist = anAxis.XDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3704 // double aYDist = anAxis.YDirection().XYZ() * ( aPnt.XYZ() - aLoc );
3705 // if ( fabs(aXDist*aXDist + aYDist*aYDist - aRad*aRad) > aToler2 )
3710 myProjector.Perform( aPnt );
3711 bool isOn = ( myProjector.IsDone() && myProjector.LowerDistance() <= myToler );
3721 ElementsOnShape::ElementsOnShape()
3723 myType(SMDSAbs_All),
3724 myToler(Precision::Confusion()),
3725 myAllNodesFlag(false)
3727 myCurShapeType = TopAbs_SHAPE;
3730 ElementsOnShape::~ElementsOnShape()
3734 void ElementsOnShape::SetMesh (const SMDS_Mesh* theMesh)
3736 myMeshModifTracer.SetMesh( theMesh );
3737 if ( myMeshModifTracer.IsMeshModified())
3738 SetShape(myShape, myType);
3741 bool ElementsOnShape::IsSatisfy (long theElementId)
3743 return myIds.Contains(theElementId);
3746 SMDSAbs_ElementType ElementsOnShape::GetType() const
3751 void ElementsOnShape::SetTolerance (const double theToler)
3753 if (myToler != theToler) {
3755 SetShape(myShape, myType);
3759 double ElementsOnShape::GetTolerance() const
3764 void ElementsOnShape::SetAllNodes (bool theAllNodes)
3766 if (myAllNodesFlag != theAllNodes) {
3767 myAllNodesFlag = theAllNodes;
3768 SetShape(myShape, myType);
3772 void ElementsOnShape::SetShape (const TopoDS_Shape& theShape,
3773 const SMDSAbs_ElementType theType)
3779 const SMDS_Mesh* myMesh = myMeshModifTracer.GetMesh();
3781 if ( !myMesh ) return;
3786 myIds.ReSize(myMesh->NbEdges() + myMesh->NbFaces() + myMesh->NbVolumes());
3789 myIds.ReSize(myMesh->NbNodes());
3792 myIds.ReSize(myMesh->NbEdges());
3795 myIds.ReSize(myMesh->NbFaces());
3797 case SMDSAbs_Volume:
3798 myIds.ReSize(myMesh->NbVolumes());
3804 myShapesMap.Clear();
3808 void ElementsOnShape::addShape (const TopoDS_Shape& theShape)
3810 if (theShape.IsNull() || myMeshModifTracer.GetMesh() == 0)
3813 if (!myShapesMap.Add(theShape)) return;
3815 myCurShapeType = theShape.ShapeType();
3816 switch (myCurShapeType)
3818 case TopAbs_COMPOUND:
3819 case TopAbs_COMPSOLID:
3823 TopoDS_Iterator anIt (theShape, Standard_True, Standard_True);
3824 for (; anIt.More(); anIt.Next()) addShape(anIt.Value());
3829 myCurSC.Load(theShape);
3835 TopoDS_Face aFace = TopoDS::Face(theShape);
3836 BRepAdaptor_Surface SA (aFace, true);
3838 u1 = SA.FirstUParameter(),
3839 u2 = SA.LastUParameter(),
3840 v1 = SA.FirstVParameter(),
3841 v2 = SA.LastVParameter();
3842 Handle(Geom_Surface) surf = BRep_Tool::Surface(aFace);
3843 myCurProjFace.Init(surf, u1,u2, v1,v2);
3850 TopoDS_Edge anEdge = TopoDS::Edge(theShape);
3851 Standard_Real u1, u2;
3852 Handle(Geom_Curve) curve = BRep_Tool::Curve(anEdge, u1, u2);
3853 myCurProjEdge.Init(curve, u1, u2);
3859 TopoDS_Vertex aV = TopoDS::Vertex(theShape);
3860 myCurPnt = BRep_Tool::Pnt(aV);
3869 void ElementsOnShape::process()
3871 const SMDS_Mesh* myMesh = myMeshModifTracer.GetMesh();
3872 if (myShape.IsNull() || myMesh == 0)
3875 if (myType == SMDSAbs_Node)
3877 SMDS_NodeIteratorPtr anIter = myMesh->nodesIterator();
3878 while (anIter->more())
3879 process(anIter->next());
3883 if (myType == SMDSAbs_Edge || myType == SMDSAbs_All)
3885 SMDS_EdgeIteratorPtr anIter = myMesh->edgesIterator();
3886 while (anIter->more())
3887 process(anIter->next());
3890 if (myType == SMDSAbs_Face || myType == SMDSAbs_All)
3892 SMDS_FaceIteratorPtr anIter = myMesh->facesIterator();
3893 while (anIter->more()) {
3894 process(anIter->next());
3898 if (myType == SMDSAbs_Volume || myType == SMDSAbs_All)
3900 SMDS_VolumeIteratorPtr anIter = myMesh->volumesIterator();
3901 while (anIter->more())
3902 process(anIter->next());
3907 void ElementsOnShape::process (const SMDS_MeshElement* theElemPtr)
3909 if (myShape.IsNull())
3912 SMDS_ElemIteratorPtr aNodeItr = theElemPtr->nodesIterator();
3913 bool isSatisfy = myAllNodesFlag;
3915 gp_XYZ centerXYZ (0, 0, 0);
3917 while (aNodeItr->more() && (isSatisfy == myAllNodesFlag))
3919 SMDS_MeshNode* aNode = (SMDS_MeshNode*)aNodeItr->next();
3920 gp_Pnt aPnt (aNode->X(), aNode->Y(), aNode->Z());
3921 centerXYZ += aPnt.XYZ();
3923 switch (myCurShapeType)
3927 myCurSC.Perform(aPnt, myToler);
3928 isSatisfy = (myCurSC.State() == TopAbs_IN || myCurSC.State() == TopAbs_ON);
3933 myCurProjFace.Perform(aPnt);
3934 isSatisfy = (myCurProjFace.IsDone() && myCurProjFace.LowerDistance() <= myToler);
3937 // check relatively the face
3938 Quantity_Parameter u, v;
3939 myCurProjFace.LowerDistanceParameters(u, v);
3940 gp_Pnt2d aProjPnt (u, v);
3941 BRepClass_FaceClassifier aClsf (myCurFace, aProjPnt, myToler);
3942 isSatisfy = (aClsf.State() == TopAbs_IN || aClsf.State() == TopAbs_ON);
3948 myCurProjEdge.Perform(aPnt);
3949 isSatisfy = (myCurProjEdge.NbPoints() > 0 && myCurProjEdge.LowerDistance() <= myToler);
3954 isSatisfy = (aPnt.Distance(myCurPnt) <= myToler);
3964 if (isSatisfy && myCurShapeType == TopAbs_SOLID) { // Check the center point for volumes MantisBug 0020168
3965 centerXYZ /= theElemPtr->NbNodes();
3966 gp_Pnt aCenterPnt (centerXYZ);
3967 myCurSC.Perform(aCenterPnt, myToler);
3968 if ( !(myCurSC.State() == TopAbs_IN || myCurSC.State() == TopAbs_ON))
3973 myIds.Add(theElemPtr->GetID());
3976 TSequenceOfXYZ::TSequenceOfXYZ()
3979 TSequenceOfXYZ::TSequenceOfXYZ(size_type n) : myArray(n)
3982 TSequenceOfXYZ::TSequenceOfXYZ(size_type n, const gp_XYZ& t) : myArray(n,t)
3985 TSequenceOfXYZ::TSequenceOfXYZ(const TSequenceOfXYZ& theSequenceOfXYZ) : myArray(theSequenceOfXYZ.myArray)
3988 template <class InputIterator>
3989 TSequenceOfXYZ::TSequenceOfXYZ(InputIterator theBegin, InputIterator theEnd): myArray(theBegin,theEnd)
3992 TSequenceOfXYZ::~TSequenceOfXYZ()
3995 TSequenceOfXYZ& TSequenceOfXYZ::operator=(const TSequenceOfXYZ& theSequenceOfXYZ)
3997 myArray = theSequenceOfXYZ.myArray;
4001 gp_XYZ& TSequenceOfXYZ::operator()(size_type n)
4003 return myArray[n-1];
4006 const gp_XYZ& TSequenceOfXYZ::operator()(size_type n) const
4008 return myArray[n-1];
4011 void TSequenceOfXYZ::clear()
4016 void TSequenceOfXYZ::reserve(size_type n)
4021 void TSequenceOfXYZ::push_back(const gp_XYZ& v)
4023 myArray.push_back(v);
4026 TSequenceOfXYZ::size_type TSequenceOfXYZ::size() const
4028 return myArray.size();
4031 TMeshModifTracer::TMeshModifTracer():
4032 myMeshModifTime(0), myMesh(0)
4035 void TMeshModifTracer::SetMesh( const SMDS_Mesh* theMesh )
4037 if ( theMesh != myMesh )
4038 myMeshModifTime = 0;
4041 bool TMeshModifTracer::IsMeshModified()
4043 bool modified = false;
4046 modified = ( myMeshModifTime != myMesh->GetMTime() );
4047 myMeshModifTime = myMesh->GetMTime();