Salome HOME
typo in comment
[modules/smesh.git] / src / StdMeshers / StdMeshers_Cartesian_3D.cxx
1 // Copyright (C) 2007-2011  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
5 //
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.
10 //
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.
15 //
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
19 //
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22 //  File   : StdMeshers_Cartesian_3D.cxx
23 //  Module : SMESH
24 //
25 #include "StdMeshers_Cartesian_3D.hxx"
26
27 #include "SMDS_MeshNode.hxx"
28 #include "SMESH_Block.hxx"
29 #include "SMESH_Comment.hxx"
30 #include "SMESH_Mesh.hxx"
31 #include "SMESH_MesherHelper.hxx"
32 #include "SMESH_subMesh.hxx"
33 #include "SMESH_subMeshEventListener.hxx"
34 #include "StdMeshers_CartesianParameters3D.hxx"
35
36 #include "utilities.h"
37 #include "Utils_ExceptHandlers.hxx"
38
39 #include <BRepAdaptor_Surface.hxx>
40 #include <BRepBndLib.hxx>
41 #include <BRepTools.hxx>
42 #include <BRep_Tool.hxx>
43 #include <Bnd_Box.hxx>
44 #include <ElSLib.hxx>
45 #include <IntAna_IntConicQuad.hxx>
46 #include <IntAna_IntLinTorus.hxx>
47 #include <IntAna_Quadric.hxx>
48 #include <IntCurveSurface_TransitionOnCurve.hxx>
49 #include <IntCurvesFace_Intersector.hxx>
50 #include <Poly_Triangulation.hxx>
51 #include <Precision.hxx>
52 #include <TopExp.hxx>
53 #include <TopExp_Explorer.hxx>
54 #include <TopTools_MapIteratorOfMapOfShape.hxx>
55 #include <TopTools_MapOfShape.hxx>
56 #include <TopoDS.hxx>
57 #include <TopoDS_Face.hxx>
58 #include <gp_Cone.hxx>
59 #include <gp_Cylinder.hxx>
60 #include <gp_Lin.hxx>
61 #include <gp_Pln.hxx>
62 #include <gp_Pnt2d.hxx>
63 #include <gp_Sphere.hxx>
64 #include <gp_Torus.hxx>
65
66 using namespace std;
67
68 #define _MY_DEBUG_
69
70 //=============================================================================
71 /*!
72  * Constructor
73  */
74 //=============================================================================
75
76 StdMeshers_Cartesian_3D::StdMeshers_Cartesian_3D(int hypId, int studyId, SMESH_Gen * gen)
77   :SMESH_3D_Algo(hypId, studyId, gen)
78 {
79   _name = "Cartesian_3D";
80   _shapeType = (1 << TopAbs_SOLID);       // 1 bit /shape type
81   _compatibleHypothesis.push_back("CartesianParameters3D");
82
83   _onlyUnaryInput = false;         // to mesh all SOLIDs at once
84   _requireDescretBoundary = false; // 2D mesh not needed
85   _supportSubmeshes = false;       // do not use any existing mesh
86 }
87
88 //=============================================================================
89 /*!
90  * Check presence of a hypothesis
91  */
92 //=============================================================================
93
94 bool StdMeshers_Cartesian_3D::CheckHypothesis (SMESH_Mesh&          aMesh,
95                                                const TopoDS_Shape&  aShape,
96                                                Hypothesis_Status&   aStatus)
97 {
98   aStatus = SMESH_Hypothesis::HYP_MISSING;
99
100   const list<const SMESHDS_Hypothesis*>& hyps = GetUsedHypothesis(aMesh, aShape);
101   list <const SMESHDS_Hypothesis* >::const_iterator h = hyps.begin();
102   if ( h == hyps.end())
103   {
104     return false;
105   }
106
107   for ( ; h != hyps.end(); ++h )
108   {
109     if (( _hyp = dynamic_cast<const StdMeshers_CartesianParameters3D*>( *h )))
110     {
111       aStatus = _hyp->IsDefined() ? HYP_OK : HYP_BAD_PARAMETER;
112       break;
113     }
114   }
115
116   return aStatus == HYP_OK;
117 }
118
119 namespace
120 {
121   //=============================================================================
122   // Definitions of internal utils
123   // --------------------------------------------------------------------------
124   enum Transition {
125     Trans_TANGENT = IntCurveSurface_Tangent,
126     Trans_IN      = IntCurveSurface_In,
127     Trans_OUT     = IntCurveSurface_Out,
128     Trans_APEX
129   };
130   /*!
131    * \brief Data of intersection between a GridLine and a TopoDS_Face
132    */
133   struct IntersectionPoint
134   {
135     double                       _paramOnLine;
136     mutable Transition           _transition;
137     mutable const SMDS_MeshNode* _node;
138
139     IntersectionPoint(): _node(0) {}
140     bool operator< ( const IntersectionPoint& o ) const { return _paramOnLine < o._paramOnLine; }
141   };
142   // --------------------------------------------------------------------------
143   /*!
144    * \brief A line of the grid and its intersections with 2D geometry
145    */
146   struct GridLine
147   {
148     gp_Lin _line;
149     double _length; // line length
150     multiset< IntersectionPoint > _intPoints;
151
152     void RemoveExcessIntPoints( const double tol );
153     bool GetIsOutBefore( multiset< IntersectionPoint >::iterator ip, bool prevIsOut );
154   };
155   // --------------------------------------------------------------------------
156   /*!
157    * \brief Iterator on the grid lines in one direction
158    */
159   struct LineIndexer
160   {
161     size_t _size  [3];
162     size_t _curInd[3];
163     size_t _iVar1, _iVar2, _iConst;
164     string _name1, _name2, _nameConst;
165     LineIndexer() {}
166     LineIndexer( size_t sz1, size_t sz2, size_t sz3,
167                  size_t iv1, size_t iv2, size_t iConst,
168                  const string& nv1, const string& nv2, const string& nConst )
169     {
170       _size[0] = sz1; _size[1] = sz2; _size[2] = sz3;
171       _curInd[0] = _curInd[1] = _curInd[2] = 0;
172       _iVar1 = iv1; _iVar2 = iv2; _iConst = iConst; 
173       _name1 = nv1; _name2 = nv2; _nameConst = nConst;
174     }
175
176     size_t I() { return _curInd[0]; }
177     size_t J() { return _curInd[1]; }
178     size_t K() { return _curInd[2]; }
179     void SetIJK( size_t i, size_t j, size_t k )
180     {
181       _curInd[0] = i; _curInd[1] = j; _curInd[2] = k;
182     }
183     void operator++()
184     {
185       if ( ++_curInd[_iVar1] == _size[_iVar1] )
186         _curInd[_iVar1] = 0, ++_curInd[_iVar2];
187     }
188     bool More() const { return _curInd[_iVar2] < _size[_iVar2]; }
189     size_t LineIndex () const { return _curInd[_iVar1] + _curInd[_iVar2]* _size[_iVar1]; }
190     size_t LineIndex10 () const { return (_curInd[_iVar1] + 1 ) + _curInd[_iVar2]* _size[_iVar1]; }
191     size_t LineIndex01 () const { return _curInd[_iVar1] + (_curInd[_iVar2] + 1 )* _size[_iVar1]; }
192     size_t LineIndex11 () const { return (_curInd[_iVar1] + 1 ) + (_curInd[_iVar2] + 1 )* _size[_iVar1]; }
193     void SetIndexOnLine (size_t i)  { _curInd[ _iConst ] = i; }
194     size_t NbLines() const { return _size[_iVar1] * _size[_iVar2]; }
195   };
196   // --------------------------------------------------------------------------
197   /*!
198    * \brief Container of GridLine's
199    */
200   struct Grid
201   {
202     vector< double >   _coords[3]; // coordinates of grid nodes
203     vector< GridLine > _lines [3]; // in 3 directions
204     double             _tol, _minCellSize;
205
206     vector< const SMDS_MeshNode* > _nodes; // mesh nodes at grid nodes
207
208     size_t NodeIndex( size_t i, size_t j, size_t k ) const
209     {
210       return i + j*_coords[0].size() + k*_coords[0].size()*_coords[1].size();
211     }
212     size_t NodeIndexDX() const { return 1; }
213     size_t NodeIndexDY() const { return _coords[0].size(); }
214     size_t NodeIndexDZ() const { return _coords[0].size() * _coords[1].size(); }
215
216     LineIndexer GetLineIndexer(size_t iDir) const;
217
218     void SetCoordinates(const vector<double>& xCoords,
219                         const vector<double>& yCoords,
220                         const vector<double>& zCoords,
221                         const TopoDS_Shape&   shape );
222     void ComputeNodes(SMESH_MesherHelper& helper);
223   };
224   // --------------------------------------------------------------------------
225   /*!
226    * \brief Intersector of TopoDS_Face with all GridLine's
227    */
228   struct FaceGridIntersector
229   {
230     TopoDS_Face _face;
231     Grid*       _grid;
232     Bnd_Box     _bndBox;
233     IntCurvesFace_Intersector* _surfaceInt;
234     vector< std::pair< GridLine*, IntersectionPoint > > _intersections;
235
236     FaceGridIntersector(): _grid(0), _surfaceInt(0) {}
237     void Intersect();
238     bool IsInGrid(const Bnd_Box& gridBox);
239
240     void StoreIntersections()
241     {
242       for ( size_t i = 0; i < _intersections.size(); ++i )
243         _intersections[i].first->_intPoints.insert( _intersections[i].second );
244     }
245     const Bnd_Box& GetFaceBndBox()
246     {
247       GetCurveFaceIntersector();
248       return _bndBox;
249     }
250     IntCurvesFace_Intersector* GetCurveFaceIntersector()
251     {
252       if ( !_surfaceInt )
253       {
254         _surfaceInt = new IntCurvesFace_Intersector( _face, Precision::PConfusion() );
255         _bndBox     = _surfaceInt->Bounding();
256         if ( _bndBox.IsVoid() )
257           BRepBndLib::Add (_face, _bndBox);
258       }
259       return _surfaceInt;
260     }
261   };
262   // --------------------------------------------------------------------------
263   /*!
264    * \brief Intersector of a surface with a GridLine
265    */
266   struct FaceLineIntersector
267   {
268     double      _tol;
269     double      _u, _v, _w; // params on the face and the line
270     Transition  _transition; // transition of at intersection (see IntCurveSurface.cdl)
271     Transition  _transIn, _transOut; // IN and OUT transitions depending of face orientation
272
273     gp_Pln      _plane;
274     gp_Cylinder _cylinder;
275     gp_Cone     _cone;
276     gp_Sphere   _sphere;
277     gp_Torus    _torus;
278     IntCurvesFace_Intersector* _surfaceInt;
279
280     vector< IntersectionPoint > _intPoints;
281
282     void IntersectWithPlane   (const GridLine& gridLine);
283     void IntersectWithCylinder(const GridLine& gridLine);
284     void IntersectWithCone    (const GridLine& gridLine);
285     void IntersectWithSphere  (const GridLine& gridLine);
286     void IntersectWithTorus   (const GridLine& gridLine);
287     void IntersectWithSurface (const GridLine& gridLine);
288
289     void addIntPoint(const bool toClassify=true);
290     bool isParamOnLineOK( const double linLength )
291     {
292       return -_tol < _w && _w < linLength + _tol;
293     }
294     FaceLineIntersector():_surfaceInt(0) {}
295     ~FaceLineIntersector() { if (_surfaceInt ) delete _surfaceInt; _surfaceInt = 0; }
296   };
297   // --------------------------------------------------------------------------
298   /*!
299    * \brief Class representing topology of the hexahedron and creating a mesh
300    *        volume basing on analysis of hexahedron intersection with geometry
301    */
302   class Hexahedron
303   {
304     // --------------------------------------------------------------------------------
305     struct _Face;
306     struct _Link;
307     // --------------------------------------------------------------------------------
308     struct _Node //!< node either at a hexahedron corner or at GridLine intersection
309     {
310       const SMDS_MeshNode*     _node; // mesh node at hexahedron corner
311       const IntersectionPoint* _intPoint;
312
313       _Node(const SMDS_MeshNode* n=0, const IntersectionPoint* ip=0):_node(n), _intPoint(ip) {} 
314       const SMDS_MeshNode* Node() const { return _intPoint ? _intPoint->_node : _node; }
315       bool IsCorner() const { return _node; }
316     };
317     // --------------------------------------------------------------------------------
318     struct _Link // link connecting two _Node's
319     {
320       _Node* _nodes[2];
321       vector< _Node>  _intNodes; // _Node's at GridLine intersections
322       vector< _Link > _splits;
323       vector< _Face*> _faces;
324     };
325     // --------------------------------------------------------------------------------
326     struct _OrientedLink
327     {
328       _Link* _link;
329       bool   _reverse;
330       _OrientedLink( _Link* link=0, bool reverse=false ): _link(link), _reverse(reverse) {}
331       void Reverse() { _reverse = !_reverse; }
332       int NbResultLinks() const { return _link->_splits.size(); }
333       _OrientedLink ResultLink(int i) const
334       {
335         return _OrientedLink(&_link->_splits[_reverse ? NbResultLinks()-i-1 : i],_reverse);
336       }
337       _Node* FirstNode() const { return _link->_nodes[ _reverse ]; }
338       _Node* LastNode() const { return _link->_nodes[ !_reverse ]; }
339       // int NbNodes() const { return 2 + _link->_intNodes.size(); }
340       // _Node* GetNode(const int i)
341       // {
342       //   return ( 0 < i && i < NbNodes()-1 ) ? _link->_intNodes[i-1] : ( _link->_nodes[bool(i)]);
343       // }
344     };
345     // --------------------------------------------------------------------------------
346     struct _Face
347     {
348       vector< _OrientedLink > _links;
349       vector< _Link >         _polyLinks; // links added to close a polygonal face
350     };
351     // --------------------------------------------------------------------------------
352     int   _nodeShift[8];
353     _Node _hexNodes[8];
354     _Link _hexLinks[12];
355     _Face _hexQuads[6];
356
357     vector< _Face > _polygons;
358
359     Grid* _grid;
360     LineIndexer _lineInd[3];
361
362     double _sizeThreshold, _sideLength[3];
363
364     int _nbCornerNodes, _nbIntNodes;
365
366   public:
367     Hexahedron(const double sizeThreshold, Grid* grid);
368     void Init( size_t i, size_t j, size_t k );
369     int MakeElements(SMESH_MesherHelper& helper);
370   private:
371     bool checkPolyhedronSize() const;
372     bool addHexa (SMESH_MesherHelper& helper);
373     bool addTetra(SMESH_MesherHelper& helper);
374     bool addPenta(SMESH_MesherHelper& helper);
375     bool addPyra (SMESH_MesherHelper& helper);
376   };
377  
378   // --------------------------------------------------------------------------
379   /*!
380    * \brief Structure intersecting certain nb of faces with GridLine's in one thread
381    */
382 #ifdef WITH_TBB
383   struct ParallelIntersector
384   {
385     vector< FaceGridIntersector >& _faceVec;
386     ParallelIntersector( vector< FaceGridIntersector >& faceVec): _faceVec(faceVec){}
387     void operator() ( const tbb::blocked_range<size_t>& r ) const
388     {
389       for ( size_t i = r.begin(); i != r.end(); ++i )
390         _faceVec[i]->Intersect();
391     }
392   };
393 #endif
394   //=============================================================================
395   // Implementation of internal utils
396   //=============================================================================
397   /*
398    * Remove coincident intersection points
399    */
400   void GridLine::RemoveExcessIntPoints( const double tol )
401   {
402     if ( _intPoints.size() < 2 ) return;
403
404     set< Transition > tranSet;
405     multiset< IntersectionPoint >::iterator ip2 = _intPoints.begin(), ip1 = ip2++;
406     for ( ; ip2 != _intPoints.end(); ip1 = ip2++ )
407     {
408       tranSet.clear();
409       while ( ip2->_paramOnLine - ip1->_paramOnLine <= tol  && ip2 != _intPoints.end())
410       {
411         tranSet.insert( ip1->_transition );
412         tranSet.insert( ip2->_transition );
413         _intPoints.erase( ip1 );
414         ip1 = ip2++;
415       }
416       if ( tranSet.size() > 1 ) // points with different transition coincide
417       {
418         bool isIN  = tranSet.count( Trans_IN );
419         bool isOUT = tranSet.count( Trans_OUT );
420         if ( isIN && isOUT )
421           (*ip1)._transition = Trans_TANGENT;
422         else
423           (*ip1)._transition = isIN ? Trans_IN : Trans_OUT;
424       }
425     }
426   }
427   //================================================================================
428   /*
429    * Return "is OUT" state for nodes before the given intersention point
430    */
431   bool GridLine::GetIsOutBefore( multiset< IntersectionPoint >::iterator ip, bool prevIsOut )
432   {
433     if ( ip->_transition == Trans_IN )
434       return true;
435     if ( ip->_transition == Trans_OUT )
436       return false;
437     if ( ip->_transition == Trans_APEX )
438     {
439       // singularity point (apex of a cone)
440       if ( _intPoints.size() == 1 || ip == _intPoints.begin() )
441         return true;
442       multiset< IntersectionPoint >::iterator ipBef = ip, ipAft = ++ip;
443       if ( ipAft == _intPoints.end() )
444         return false;
445       --ipBef;
446       if ( ipBef->_transition != ipAft->_transition )
447         return ( ipBef->_transition == Trans_OUT );
448       return ( ipBef->_transition != Trans_OUT );
449     }
450     return prevIsOut; // _transition == Trans_TANGENT
451   }
452   //================================================================================
453   /*
454    * Return an iterator on GridLine's in a given direction
455    */
456   LineIndexer Grid::GetLineIndexer(size_t iDir) const
457   {
458     const size_t indices[] = { 1,2,0, 0,2,1, 0,1,2 };
459     const string s[] = { "X", "Y", "Z" };
460     LineIndexer li( _coords[0].size(),  _coords[1].size(),    _coords[2].size(),
461                     indices[iDir*3],    indices[iDir*3+1],    indices[iDir*3+2],
462                     s[indices[iDir*3]], s[indices[iDir*3+1]], s[indices[iDir*3+2]]);
463     return li;
464   }
465   //=============================================================================
466   /*
467    * Creates GridLine's of the grid
468    */
469   void Grid::SetCoordinates(const vector<double>& xCoords,
470                             const vector<double>& yCoords,
471                             const vector<double>& zCoords,
472                             const TopoDS_Shape&   shape)
473   {
474     _coords[0] = xCoords;
475     _coords[1] = yCoords;
476     _coords[2] = zCoords;
477
478     // compute tolerance
479     _minCellSize = Precision::Infinite();
480     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
481     {
482       for ( size_t i = 1; i < _coords[ iDir ].size(); ++i )
483       {
484         double cellLen = _coords[ iDir ][ i ] - _coords[ iDir ][ i-1 ];
485         if ( cellLen < _minCellSize )
486           _minCellSize = cellLen;
487       }
488     }
489     if ( _minCellSize < Precision::Confusion() )
490       throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
491                                 SMESH_Comment("Too small cell size: ") << _tol );
492     _tol = _minCellSize / 1000.;
493
494     // attune grid extremities to shape bounding box computed by vertices
495     Bnd_Box shapeBox;
496     for ( TopExp_Explorer vExp( shape, TopAbs_VERTEX ); vExp.More(); vExp.Next() )
497       shapeBox.Add( BRep_Tool::Pnt( TopoDS::Vertex( vExp.Current() )));
498     
499     double sP[6]; // aXmin, aYmin, aZmin, aXmax, aYmax, aZmax
500     shapeBox.Get(sP[0],sP[1],sP[2],sP[3],sP[4],sP[5]);
501     double* cP[6] = { &_coords[0].front(), &_coords[1].front(), &_coords[2].front(),
502                       &_coords[0].back(),  &_coords[1].back(),  &_coords[2].back() };
503     for ( int i = 0; i < 6; ++i )
504       if ( fabs( sP[i] - *cP[i] ) < _tol )
505         *cP[i] = sP[i] + _tol/1000. * ( i < 3 ? +1 : -1 );
506
507     // create lines
508     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
509     {
510       LineIndexer li = GetLineIndexer( iDir );
511       _lines[iDir].resize( li.NbLines() );
512       double len = _coords[ iDir ].back() - _coords[iDir].front();
513       gp_Vec dir( iDir==0, iDir==1, iDir==2 );
514       for ( ; li.More(); ++li )
515       {
516         GridLine& gl = _lines[iDir][ li.LineIndex() ];
517         gl._line.SetLocation(gp_Pnt(_coords[0][li.I()], _coords[1][li.J()], _coords[2][li.K()])); 
518         gl._line.SetDirection( dir );
519         gl._length = len;
520       }
521     }
522   }
523   //================================================================================
524   /*
525    * Creates all nodes
526    */
527   void Grid::ComputeNodes(SMESH_MesherHelper& helper)
528   {
529     // state of each node of the grid relative to the geomerty
530     vector< bool > isNodeOut( _coords[0].size() * _coords[1].size() * _coords[2].size(), false );
531     _nodes.resize( isNodeOut.size(), 0 );
532
533     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
534     {
535       LineIndexer li = GetLineIndexer( iDir );
536
537       // find out a shift of node index while walking along a GridLine in this direction
538       li.SetIndexOnLine( 0 );
539       size_t nIndex0 = NodeIndex( li.I(), li.J(), li.K() );
540       li.SetIndexOnLine( 1 );
541       const size_t nShift = NodeIndex( li.I(), li.J(), li.K() ) - nIndex0;
542       
543       const vector<double> & coords = _coords[ iDir ];
544       for ( ; li.More(); ++li ) // loop on lines in iDir
545       {
546         li.SetIndexOnLine( 0 );
547         nIndex0 = NodeIndex( li.I(), li.J(), li.K() );
548
549         GridLine& line = _lines[ iDir ][ li.LineIndex() ];
550         line.RemoveExcessIntPoints( _tol );
551         multiset< IntersectionPoint >& intPnts = _lines[ iDir ][ li.LineIndex() ]._intPoints;
552         multiset< IntersectionPoint >::iterator ip = intPnts.begin();
553
554         bool isOut = true;
555         const double* nodeCoord = & coords[0], *coord0 = nodeCoord, *coordEnd = coord0 + coords.size();
556         double nodeParam = 0;
557         for ( ; ip != intPnts.end(); ++ip )
558         {
559           // set OUT state or just skip IN nodes before ip
560           if ( nodeParam < ip->_paramOnLine - _tol )
561           {
562             isOut = line.GetIsOutBefore( ip, isOut );
563
564             while ( nodeParam < ip->_paramOnLine - _tol )
565             {
566               if ( isOut )
567                 isNodeOut[ nIndex0 + nShift * ( nodeCoord-coord0 ) ] = isOut;
568               if ( ++nodeCoord <  coordEnd )
569                 nodeParam = *nodeCoord - *coord0;
570               else
571                 break;
572             }
573             if ( nodeCoord == coordEnd ) break;
574           }
575           // create a mesh node on a GridLine at ip if it does not coincide with a grid node
576           if ( nodeParam > ip->_paramOnLine + _tol )
577           {
578             li.SetIndexOnLine( 0 );
579             double xyz[3] = { _coords[0][ li.I() ], _coords[1][ li.J() ], _coords[2][ li.K() ]};
580             xyz[ li._iConst ] += ip->_paramOnLine;
581             ip->_node = helper.AddNode( xyz[0], xyz[1], xyz[2] );
582           }
583           // create a mesh node at ip concident with a grid node
584           else
585           {
586             int nodeIndex = nIndex0 + nShift * ( nodeCoord-coord0 );
587             if ( ! _nodes[ nodeIndex ] )
588             {
589               li.SetIndexOnLine( nodeCoord-coord0 );
590               double xyz[3] = { _coords[0][ li.I() ], _coords[1][ li.J() ], _coords[2][ li.K() ]};
591               _nodes[ nodeIndex ] = helper.AddNode( xyz[0], xyz[1], xyz[2] );
592             }
593             if ( ++nodeCoord < coordEnd )
594               nodeParam = *nodeCoord - *coord0;
595           }
596         }
597         // set OUT state to nodes after the last ip
598         for ( ; nodeCoord < coordEnd; ++nodeCoord )
599           isNodeOut[ nIndex0 + nShift * ( nodeCoord-coord0 ) ] = true;
600       }
601     }
602
603     // Create mesh nodes at !OUT nodes of the grid
604
605     for ( size_t z = 0; z < _coords[2].size(); ++z )
606       for ( size_t y = 0; y < _coords[1].size(); ++y )
607         for ( size_t x = 0; x < _coords[0].size(); ++x )
608         {
609           size_t nodeIndex = NodeIndex( x, y, z );
610           if ( !isNodeOut[ nodeIndex ] && !_nodes[ nodeIndex] )
611             _nodes[ nodeIndex ] = helper.AddNode( _coords[0][x], _coords[1][y], _coords[2][z] );
612         }
613
614 #ifdef _MY_DEBUG_
615     // check validity of transitions
616     const char* trName[] = { "TANGENT", "IN", "OUT", "APEX" };
617     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
618     {
619       LineIndexer li = GetLineIndexer( iDir );
620       for ( ; li.More(); ++li )
621       {
622         multiset< IntersectionPoint >& intPnts = _lines[ iDir ][ li.LineIndex() ]._intPoints;
623         if ( intPnts.empty() ) continue;
624         if ( intPnts.size() == 1 )
625         {
626           if ( intPnts.begin()->_transition != Trans_TANGENT &&
627                intPnts.begin()->_transition != Trans_APEX )
628           throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
629                                     SMESH_Comment("Wrong SOLE transition of GridLine (")
630                                     << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
631                                     << ") along " << li._nameConst
632                                     << ": " << trName[ intPnts.begin()->_transition] );
633         }
634         else
635         {
636           if ( intPnts.begin()->_transition == Trans_OUT )
637             throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
638                                       SMESH_Comment("Wrong START transition of GridLine (")
639                                       << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
640                                       << ") along " << li._nameConst
641                                       << ": " << trName[ intPnts.begin()->_transition ]);
642           if ( intPnts.rbegin()->_transition == Trans_IN )
643             throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
644                                       SMESH_Comment("Wrong END transition of GridLine (")
645                                       << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
646                                       << ") along " << li._nameConst
647                                     << ": " << trName[ intPnts.rbegin()->_transition ]);
648         }
649       }
650     }
651 #endif
652   }
653
654   //=============================================================================
655   /*
656    * Checks if the face is encosed by the grid
657    */
658   bool FaceGridIntersector::IsInGrid(const Bnd_Box& gridBox)
659   {
660     double x0,y0,z0, x1,y1,z1;
661     const Bnd_Box& faceBox = GetFaceBndBox();
662     faceBox.Get(x0,y0,z0, x1,y1,z1);
663
664     if ( !gridBox.IsOut( gp_Pnt( x0,y0,z0 )) &&
665          !gridBox.IsOut( gp_Pnt( x1,y1,z1 )))
666       return true;
667
668     double X0,Y0,Z0, X1,Y1,Z1;
669     gridBox.Get(X0,Y0,Z0, X1,Y1,Z1);
670     double faceP[6] = { x0,y0,z0, x1,y1,z1 };
671     double gridP[6] = { X0,Y0,Z0, X1,Y1,Z1 };
672     gp_Dir axes[3]  = { gp::DX(), gp::DY(), gp::DZ() };
673     for ( int iDir = 0; iDir < 6; ++iDir )
674     {
675       if ( iDir < 3  && gridP[ iDir ] <= faceP[ iDir ] ) continue;
676       if ( iDir >= 3 && gridP[ iDir ] >= faceP[ iDir ] ) continue;
677
678       // check if the face intersects a side of a gridBox
679
680       gp_Pnt p = iDir < 3 ? gp_Pnt( X0,Y0,Z0 ) : gp_Pnt( X1,Y1,Z1 );
681       gp_Ax1 norm( p, axes[ iDir % 3 ] );
682       if ( iDir < 3 ) norm.Reverse();
683
684       gp_XYZ O = norm.Location().XYZ(), N = norm.Direction().XYZ();
685
686       TopLoc_Location loc = _face.Location();
687       Handle(Poly_Triangulation) aPoly = BRep_Tool::Triangulation(_face,loc);
688       if ( !aPoly.IsNull() )
689       {
690         if ( !loc.IsIdentity() )
691         {
692           norm.Transform( loc.Transformation().Inverted() );
693           O = norm.Location().XYZ(), N = norm.Direction().XYZ();
694         }
695         const double deflection = aPoly->Deflection();
696
697         const TColgp_Array1OfPnt& nodes = aPoly->Nodes();
698         for ( int i = nodes.Lower(); i <= nodes.Upper(); ++i )
699           if (( nodes( i ).XYZ() - O ) * N > _grid->_tol + deflection )
700             return false;
701       }
702       else
703       {
704         BRepAdaptor_Surface surf( _face );
705         double u0, u1, v0, v1, du, dv, u, v;
706         BRepTools::UVBounds( _face, u0, u1, v0, v1);
707         if ( surf.GetType() == GeomAbs_Plane ) {
708           du = u1 - u0, dv = v1 - v0;
709         }
710         else {
711           du = surf.UResolution( _grid->_minCellSize / 10. );
712           dv = surf.VResolution( _grid->_minCellSize / 10. );
713         }
714         for ( u = u0, v = v0; u <= u1 && v <= v1; u += du, v += dv )
715         {
716           gp_Pnt p = surf.Value( u, v );
717           if (( p.XYZ() - O ) * N > _grid->_tol )
718           {
719             TopAbs_State state = GetCurveFaceIntersector()->ClassifyUVPoint(gp_Pnt2d( u, v ));
720             if ( state == TopAbs_IN || state == TopAbs_ON )
721               return false;
722           }
723         }
724       }
725     }
726     return true;
727   }
728   //=============================================================================
729   /*
730    * Intersects TopoDS_Face with all GridLine's
731    */
732   void FaceGridIntersector::Intersect()
733   {
734     FaceLineIntersector intersector;
735     intersector._surfaceInt = GetCurveFaceIntersector();
736     intersector._tol        = _grid->_tol;
737     intersector._transOut   = _face.Orientation() == TopAbs_REVERSED ? Trans_IN : Trans_OUT;
738     intersector._transIn    = _face.Orientation() == TopAbs_REVERSED ? Trans_OUT : Trans_IN;
739
740     typedef void (FaceLineIntersector::* PIntFun )(const GridLine& gridLine);
741     PIntFun interFunction;
742
743     BRepAdaptor_Surface surf( _face );
744     switch ( surf.GetType() ) {
745     case GeomAbs_Plane:
746       intersector._plane = surf.Plane();
747       interFunction = &FaceLineIntersector::IntersectWithPlane;
748       break;
749     case GeomAbs_Cylinder:
750       intersector._cylinder = surf.Cylinder();
751       interFunction = &FaceLineIntersector::IntersectWithCylinder;
752       break;
753     case GeomAbs_Cone:
754       intersector._cone = surf.Cone();
755       interFunction = &FaceLineIntersector::IntersectWithCone;
756       break;
757     case GeomAbs_Sphere:
758       intersector._sphere = surf.Sphere();
759       interFunction = &FaceLineIntersector::IntersectWithSphere;
760       break;
761     case GeomAbs_Torus:
762       intersector._torus = surf.Torus();
763       interFunction = &FaceLineIntersector::IntersectWithTorus;
764       break;
765     default:
766       interFunction = &FaceLineIntersector::IntersectWithSurface;
767     }
768
769     _intersections.clear();
770     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
771     {
772       if ( surf.GetType() == GeomAbs_Plane )
773       {
774         // check if all lines in this direction are parallel to a plane
775         if ( intersector._plane.Axis().IsNormal( _grid->_lines[iDir][0]._line.Position(),
776                                                  Precision::Angular()))
777           continue;
778         // find out a transition, that is the same for all lines of a direction
779         gp_Dir plnNorm = intersector._plane.Axis().Direction();
780         gp_Dir lineDir = _grid->_lines[iDir][0]._line.Direction();
781         intersector._transition =
782           ( plnNorm * lineDir < 0 ) ? intersector._transIn : intersector._transOut;
783       }
784       if ( surf.GetType() == GeomAbs_Cylinder )
785       {
786         // check if all lines in this direction are parallel to a cylinder
787         if ( intersector._cylinder.Axis().IsParallel( _grid->_lines[iDir][0]._line.Position(),
788                                                       Precision::Angular()))
789           continue;
790       }
791
792       // intersect the grid lines with the face
793       for ( size_t iL = 0; iL < _grid->_lines[iDir].size(); ++iL )
794       {
795         GridLine& gridLine = _grid->_lines[iDir][iL];
796         if ( _bndBox.IsOut( gridLine._line )) continue;
797
798         intersector._intPoints.clear();
799         (intersector.*interFunction)( gridLine );
800         for ( size_t i = 0; i < intersector._intPoints.size(); ++i )
801           _intersections.push_back( make_pair( &gridLine, intersector._intPoints[i] ));
802       }
803     }
804   }
805   //================================================================================
806   /*
807    * Store an intersection if it is In or ON the face
808    */
809   void FaceLineIntersector::addIntPoint(const bool toClassify)
810   {
811     TopAbs_State state = toClassify ? _surfaceInt->ClassifyUVPoint(gp_Pnt2d( _u, _v )) : TopAbs_IN;
812     if ( state == TopAbs_IN || state == TopAbs_ON )
813     {
814       IntersectionPoint p;
815       p._paramOnLine = _w;
816       p._transition  = _transition;
817       _intPoints.push_back( p );
818     }
819   }
820   //================================================================================
821   /*
822    * Intersect a line with a plane
823    */
824   void FaceLineIntersector::IntersectWithPlane   (const GridLine& gridLine)
825   {
826     IntAna_IntConicQuad linPlane( gridLine._line, _plane, Precision::Angular());
827     _w = linPlane.ParamOnConic(1);
828     if ( isParamOnLineOK( gridLine._length ))
829     {
830       ElSLib::Parameters(_plane, linPlane.Point(1) ,_u,_v);
831       addIntPoint();
832     }
833   }
834   //================================================================================
835   /*
836    * Intersect a line with a cylinder
837    */
838   void FaceLineIntersector::IntersectWithCylinder(const GridLine& gridLine)
839   {
840     IntAna_IntConicQuad linCylinder( gridLine._line,_cylinder);
841     if ( linCylinder.IsDone() && linCylinder.NbPoints() > 0 )
842     {
843       _w = linCylinder.ParamOnConic(1);
844       if ( linCylinder.NbPoints() == 1 )
845         _transition = Trans_TANGENT;
846       else
847         _transition = _w < linCylinder.ParamOnConic(2) ? _transIn : _transOut;
848       if ( isParamOnLineOK( gridLine._length ))
849       {
850         ElSLib::Parameters(_cylinder, linCylinder.Point(1) ,_u,_v);
851         addIntPoint();
852       }
853       if ( linCylinder.NbPoints() > 1 )
854       {
855         _w = linCylinder.ParamOnConic(2);
856         if ( isParamOnLineOK( gridLine._length ))
857         {
858           ElSLib::Parameters(_cylinder, linCylinder.Point(2) ,_u,_v);
859           _transition = ( _transition == Trans_OUT ) ? Trans_IN : Trans_OUT;
860           addIntPoint();
861         }
862       }
863     }
864   }
865   //================================================================================
866   /*
867    * Intersect a line with a cone
868    */
869   void FaceLineIntersector::IntersectWithCone (const GridLine& gridLine)
870   {
871     IntAna_IntConicQuad linCone(gridLine._line,_cone);
872     if ( !linCone.IsDone() ) return;
873     gp_Pnt P;
874     gp_Vec du, dv, norm;
875     for ( int i = 1; i <= linCone.NbPoints(); ++i )
876     {
877       _w = linCone.ParamOnConic( i );
878       if ( !isParamOnLineOK( gridLine._length )) continue;
879       ElSLib::Parameters(_cone, linCone.Point(i) ,_u,_v);
880       TopAbs_State state = _surfaceInt->ClassifyUVPoint(gp_Pnt2d( _u, _v ));
881       if ( state == TopAbs_IN || state == TopAbs_ON )
882       {
883         ElSLib::D1( _u, _v, _cone, P, du, dv );
884         norm = du ^ dv;
885         double normSize2 = norm.SquareMagnitude();
886         if ( normSize2 > Precision::Angular() * Precision::Angular() )
887         {
888           double cos = norm.XYZ() * gridLine._line.Direction().XYZ();
889           cos /= sqrt( normSize2 );
890           if ( cos < -Precision::Angular() )
891             _transition = _transIn;
892           else if ( cos > Precision::Angular() )
893             _transition = _transOut;
894           else
895             _transition = Trans_TANGENT;
896         }
897         else
898         {
899           _transition = Trans_APEX;
900         }
901         addIntPoint( /*toClassify=*/false);
902       }
903     }
904   }
905   //================================================================================
906   /*
907    * Intersect a line with a sphere
908    */
909   void FaceLineIntersector::IntersectWithSphere  (const GridLine& gridLine)
910   {
911     IntAna_IntConicQuad linSphere(gridLine._line,_sphere);
912     if ( linSphere.IsDone() && linSphere.NbPoints() > 0 )
913     {
914       _w = linSphere.ParamOnConic(1);
915       if ( linSphere.NbPoints() == 1 )
916         _transition = Trans_TANGENT;
917       else
918         _transition = _w < linSphere.ParamOnConic(2) ? _transIn : _transOut;
919       if ( isParamOnLineOK( gridLine._length ))
920       {
921         ElSLib::Parameters(_sphere, linSphere.Point(1) ,_u,_v);
922         addIntPoint();
923       }
924       if ( linSphere.NbPoints() > 1 )
925       {
926         _w = linSphere.ParamOnConic(2);
927         if ( isParamOnLineOK( gridLine._length ))
928         {
929           ElSLib::Parameters(_sphere, linSphere.Point(2) ,_u,_v);
930           _transition = ( _transition == Trans_OUT ) ? Trans_IN : Trans_OUT;
931           addIntPoint();
932         }
933       }
934     }
935   }
936   //================================================================================
937   /*
938    * Intersect a line with a torus
939    */
940   void FaceLineIntersector::IntersectWithTorus   (const GridLine& gridLine)
941   {
942     IntAna_IntLinTorus linTorus(gridLine._line,_torus);
943     if ( !linTorus.IsDone()) return;
944     gp_Pnt P;
945     gp_Vec du, dv, norm;
946     for ( int i = 1; i <= linTorus.NbPoints(); ++i )
947     {
948       _w = linTorus.ParamOnLine( i );
949       if ( !isParamOnLineOK( gridLine._length )) continue;
950       linTorus.ParamOnTorus( i, _u,_v );
951       TopAbs_State state = _surfaceInt->ClassifyUVPoint(gp_Pnt2d( _u, _v ));
952       if ( state == TopAbs_IN || state == TopAbs_ON )
953       {
954         ElSLib::D1( _u, _v, _torus, P, du, dv );
955         norm = du ^ dv;
956         double normSize = norm.Magnitude();
957         double cos = norm.XYZ() * gridLine._line.Direction().XYZ();
958         cos /= normSize;
959         if ( cos < -Precision::Angular() )
960           _transition = _transIn;
961         else if ( cos > Precision::Angular() )
962           _transition = _transOut;
963         else
964           _transition = Trans_TANGENT;
965         addIntPoint( /*toClassify=*/false);
966       }
967     }
968   }
969   //================================================================================
970   /*
971    * Intersect a line with a non-analytical surface
972    */
973   void FaceLineIntersector::IntersectWithSurface (const GridLine& gridLine)
974   {
975     _surfaceInt->Perform( gridLine._line, 0.0, gridLine._length );
976     if ( !_surfaceInt->IsDone() ) return;
977     for ( int i = 1; i <= _surfaceInt->NbPnt(); ++i )
978     {
979       _transition = Transition( _surfaceInt->Transition( i ) );
980       _w = _surfaceInt->WParameter( i );
981       addIntPoint(/*toClassify=*/false);
982     }
983   }
984
985   //================================================================================
986   /*!
987    * \brief Creates topology of the hexahedron
988    */
989   Hexahedron::Hexahedron(const double sizeThreshold, Grid* grid)
990     : _grid( grid ), _sizeThreshold(sizeThreshold)
991   {
992     _lineInd[0] = grid->GetLineIndexer( 0 );
993     _lineInd[1] = grid->GetLineIndexer( 1 );
994     _lineInd[2] = grid->GetLineIndexer( 2 );
995
996     _polygons.reserve(100); // to avoid reallocation;
997
998     //set nodes shift within grid->_nodes from the node 000 
999     size_t dx = _grid->NodeIndexDX();
1000     size_t dy = _grid->NodeIndexDY();
1001     size_t dz = _grid->NodeIndexDZ();
1002     size_t i000 = 0;
1003     size_t i100 = i000 + dx;
1004     size_t i010 = i000 + dy;
1005     size_t i110 = i010 + dx;
1006     size_t i001 = i000 + dz;
1007     size_t i101 = i100 + dz;
1008     size_t i011 = i010 + dz;
1009     size_t i111 = i110 + dz;
1010     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V000 )] = i000;
1011     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V100 )] = i100;
1012     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V010 )] = i010;
1013     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V110 )] = i110;
1014     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V001 )] = i001;
1015     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V101 )] = i101;
1016     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V011 )] = i011;
1017     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V111 )] = i111;
1018
1019     vector< int > idVec;
1020     // set nodes to links
1021     for ( int linkID = SMESH_Block::ID_Ex00; linkID <= SMESH_Block::ID_E11z; ++linkID )
1022     {
1023       SMESH_Block::GetEdgeVertexIDs( linkID, idVec );
1024       _Link& link = _hexLinks[ SMESH_Block::ShapeIndex( linkID )];
1025       link._nodes[0] = &_hexNodes[ SMESH_Block::ShapeIndex( idVec[0] )];
1026       link._nodes[1] = &_hexNodes[ SMESH_Block::ShapeIndex( idVec[1] )];
1027       link._intNodes.reserve( 10 ); // to avoid reallocation
1028       link._splits.reserve( 10 );
1029     }
1030
1031     // set links to faces
1032     int interlace[4] = { 0, 3, 1, 2 }; // to walk by links around a face: { u0, 1v, u1, 0v }
1033     for ( int faceID = SMESH_Block::ID_Fxy0; faceID <= SMESH_Block::ID_F1yz; ++faceID )
1034     {
1035       SMESH_Block::GetFaceEdgesIDs( faceID, idVec );
1036       _Face& quad = _hexQuads[ SMESH_Block::ShapeIndex( faceID )];
1037       bool revFace = ( faceID == SMESH_Block::ID_Fxy0 ||
1038                        faceID == SMESH_Block::ID_Fx1z ||
1039                        faceID == SMESH_Block::ID_F0yz );
1040       quad._links.resize(4);
1041       vector<_OrientedLink>::iterator         frwLinkIt = quad._links.begin();
1042       vector<_OrientedLink>::reverse_iterator revLinkIt = quad._links.rbegin();
1043       for ( int i = 0; i < 4; ++i )
1044       {
1045         bool revLink = revFace;
1046         if ( i > 1 ) // to reverse u1 and v0
1047           revLink = !revLink;
1048         _OrientedLink& link = revFace ? *revLinkIt++ : *frwLinkIt++;
1049         link = _OrientedLink( & _hexLinks[ SMESH_Block::ShapeIndex( idVec[interlace[i]] )],
1050                               revLink );
1051       }
1052     }
1053   }
1054   //================================================================================
1055   /*!
1056    * \brief Initializes its data by given grid cell
1057    */
1058   void Hexahedron::Init( size_t i, size_t j, size_t k )
1059   {
1060     // set nodes of grid to nodes of the hexahedron and
1061     // count nodes at hexahedron corners located IN geometry
1062     _nbCornerNodes = _nbIntNodes = 0;
1063     size_t i000 = _grid->NodeIndex( i,j,k );
1064     for ( int iN = 0; iN < 8; ++iN )
1065       _nbCornerNodes += bool(( _hexNodes[iN]._node = _grid->_nodes[ i000 + _nodeShift[iN] ]));
1066
1067     // set intersection nodes from GridLine's to hexahedron links
1068     int linkID = 0;
1069     _Link split;
1070     IntersectionPoint curIntPnt;
1071     size_t ijk[3] = { i, j, k };
1072     for ( int iDir = 0; iDir < 3; ++iDir )
1073     {
1074       _lineInd[ iDir ].SetIJK( i,j,k );
1075       size_t lineIndex[4] = {
1076         _lineInd[ iDir ].LineIndex(),
1077         _lineInd[ iDir ].LineIndex10(),
1078         _lineInd[ iDir ].LineIndex01(),
1079         _lineInd[ iDir ].LineIndex11()
1080       };
1081       const vector<double>& coords = _grid->_coords[ iDir ];
1082       double nodeParam1 = coords[ ijk[ iDir ]    ] - coords[0] + _grid->_tol;
1083       double nodeParam2 = coords[ ijk[ iDir ] + 1] - coords[0] - _grid->_tol;
1084       _sideLength[ iDir ] = nodeParam2 - nodeParam1 + 2 * _grid->_tol;      
1085       for ( int iL = 0; iL < 4; ++iL  )
1086       {
1087         GridLine& line = _grid->_lines[ iDir ][ lineIndex[ iL ]];
1088         _Link&    link = _hexLinks[ linkID++ ];
1089         link._intNodes.clear();
1090         link._splits.clear();
1091         split._nodes[ 0 ] = link._nodes[0];
1092         curIntPnt._paramOnLine = nodeParam1;
1093         multiset< IntersectionPoint >::const_iterator ip = line._intPoints.lower_bound( curIntPnt );
1094         while ( ip != line._intPoints.end() &&
1095                 ip->_paramOnLine <= nodeParam2 &&
1096                 ip->_node )
1097         {
1098           link._intNodes.push_back( _Node( 0, &(*ip) ));
1099           ++_nbIntNodes;
1100           ++ip;
1101           // create sub-links (_splits) by splitting a link with _intNodes
1102           if ( split._nodes[ 0 ]->Node() )
1103           {
1104             split._nodes[ 1 ] = &link._intNodes.back();
1105             link._splits.push_back( split );
1106           }
1107           split._nodes[ 0 ] = &link._intNodes.back();
1108         }
1109         if ( link._nodes[ 1 ]->Node() && split._nodes[ 0 ]->Node() )
1110         {
1111           split._nodes[ 1 ] = link._nodes[1];
1112           link._splits.push_back( split );
1113         }
1114       }
1115     }
1116   }
1117
1118   //================================================================================
1119   /*!
1120    * \brief Creates mesh volumes
1121    */
1122   int Hexahedron::MakeElements(SMESH_MesherHelper& helper)
1123   {
1124     int nbAddedVols = 0;
1125     if ( _nbCornerNodes == 8 && _nbIntNodes == 0 )
1126     {
1127       // order of _hexNodes is defined by enum SMESH_Block::TShapeID
1128       helper.AddVolume( _hexNodes[0].Node(), _hexNodes[2].Node(),
1129                         _hexNodes[3].Node(), _hexNodes[1].Node(),
1130                         _hexNodes[4].Node(), _hexNodes[6].Node(),
1131                         _hexNodes[7].Node(), _hexNodes[5].Node() );
1132       return 1;
1133     }
1134     if ( _nbCornerNodes + _nbIntNodes < 4 )
1135       return nbAddedVols;
1136
1137     _polygons.clear();
1138
1139     vector<const SMDS_MeshNode* > polyhedraNodes;
1140     vector<int>                   quantities;
1141
1142     // create polygons from quadrangles and get their nodes
1143
1144     vector<_Node*> nodes;
1145     nodes.reserve( _nbCornerNodes + _nbIntNodes );
1146
1147     _Link polyLink;
1148     polyLink._faces.reserve( 1 );
1149
1150     for ( int iF = 0; iF < 6; ++iF )
1151     {
1152       const _Face& quad = _hexQuads[ iF ] ;
1153
1154       _polygons.resize( _polygons.size() + 1 );
1155       _Face& polygon = _polygons.back();
1156       polygon._links.clear();
1157       polygon._polyLinks.clear(); polygon._polyLinks.reserve( 10 );
1158
1159       // add splits of a link to a polygon and collect info on nodes
1160       //int nbIn = 0, nbOut = 0, nbCorners = 0;
1161       nodes.clear();
1162       for ( int iE = 0; iE < 4; ++iE )
1163       {
1164         int nbSpits = quad._links[ iE ].NbResultLinks();
1165         for ( int iS = 0; iS < nbSpits; ++iS )
1166         {
1167           _OrientedLink split = quad._links[ iE ].ResultLink( iS );
1168           _Node* n = split.FirstNode();
1169           if ( !polygon._links.empty() )
1170           {
1171             _Node* nPrev = polygon._links.back().LastNode();
1172             if ( nPrev != n )
1173             {
1174               polyLink._nodes[0] = nPrev;
1175               polyLink._nodes[1] = n;
1176               polygon._polyLinks.push_back( polyLink );
1177               polygon._links.push_back( _OrientedLink( &polygon._polyLinks.back() ));
1178               nodes.push_back( nPrev );
1179             }
1180           }
1181           polygon._links.push_back( split );
1182           nodes.push_back( n );
1183
1184           // if ( n->IsCorner() )
1185           //   ++nbCorners;
1186           // if ( n->_intPoint )
1187           // {
1188           //   if ( n->_intPoint->_transition == Trans_IN )
1189           //     ++nbIn;
1190           //   else if ( n->_intPoint->_transition == Trans_OUT )
1191           //     ++nbOut;
1192           //   else
1193           //     ++nbIn, ++nbOut;
1194           // }
1195         }
1196       }
1197       if ( polygon._links.size() > 1 )
1198       {
1199         _Node* n1 = polygon._links.back().LastNode();
1200         _Node* n2 = polygon._links.front().FirstNode();
1201         if ( n1 != n2 )
1202         {
1203           polyLink._nodes[0] = n1;
1204           polyLink._nodes[1] = n2;
1205           polygon._polyLinks.push_back( polyLink );
1206           polygon._links.push_back( _OrientedLink( &polygon._polyLinks.back() ));
1207           nodes.push_back( n1 );
1208         }
1209         // add polygon to its links
1210         for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
1211           polygon._links[ iL ]._link->_faces.push_back( &polygon );
1212         // store polygon nodes
1213         quantities.push_back( nodes.size() );
1214         for ( size_t i = 0; i < nodes.size(); ++i )
1215           polyhedraNodes.push_back( nodes[i]->Node() );
1216       }
1217       else
1218       {
1219         _polygons.resize( _polygons.size() - 1 );
1220       }
1221     }
1222
1223     // create polygons closing holes in a polyhedron
1224
1225     // find free links
1226     vector< _OrientedLink* > freeLinks;
1227     for ( size_t iP = 0; iP < _polygons.size(); ++iP )
1228     {
1229       _Face& polygon = _polygons[ iP ];
1230       for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
1231         if ( polygon._links[ iL ]._link->_faces.size() < 2 )
1232           freeLinks.push_back( & polygon._links[ iL ]);
1233     }
1234     // make closed chains of free links
1235     int nbFreeLinks = freeLinks.size();
1236     if ( 0 < nbFreeLinks && nbFreeLinks < 3 ) return nbAddedVols;
1237     while ( nbFreeLinks > 0 )
1238     {
1239       nodes.clear();
1240       _polygons.resize( _polygons.size() + 1 );
1241       _Face& polygon = _polygons.back();
1242       polygon._links.clear();
1243
1244       // get a remaining link to start from
1245       _OrientedLink* curLink = 0;
1246       for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
1247         if (( curLink = freeLinks[ iL ] ))
1248           freeLinks[ iL ] = 0;
1249       nodes.push_back( curLink->LastNode() );
1250       polygon._links.push_back( *curLink );
1251
1252       // find all links connected to curLink
1253       _Node* curNode = 0;
1254       do
1255       {
1256         curNode = curLink->FirstNode();
1257         curLink = 0;
1258         for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
1259           if ( freeLinks[ iL ] && freeLinks[ iL ]->LastNode() == curNode )
1260           {
1261             curLink = freeLinks[ iL ];
1262             freeLinks[ iL ] = 0;
1263             nodes.push_back( curNode );
1264             polygon._links.push_back( *curLink );
1265           }
1266       } while ( curLink );
1267
1268       nbFreeLinks -= polygon._links.size();
1269
1270       if ( curNode != nodes.front() || polygon._links.size() < 3 )
1271         return nbAddedVols; // closed polygon not found -> invalid polyhedron
1272
1273       quantities.push_back( nodes.size() );
1274       for ( size_t i = 0; i < nodes.size(); ++i )
1275         polyhedraNodes.push_back( nodes[i]->Node() );
1276
1277       // add polygon to its links and reverse links
1278       for ( size_t i = 0; i < polygon._links.size(); ++i )
1279       {
1280         polygon._links[i].Reverse();
1281         polygon._links[i]._link->_faces.push_back( &polygon );
1282       }
1283
1284       //const size_t firstPoly = _polygons.size();
1285     }
1286
1287     if ( ! checkPolyhedronSize() )
1288       return nbAddedVols;
1289
1290     // create a classic cell if possible
1291     const int nbNodes = _nbCornerNodes + _nbIntNodes;
1292     if (      nbNodes == 8 && _polygons.size() == 6 && addHexa ( helper ))
1293       ++nbAddedVols;
1294     else if ( nbNodes == 4 && _polygons.size() == 4 && addTetra( helper ))
1295       ++nbAddedVols;
1296     else if ( nbNodes == 6 && _polygons.size() == 5 && addPenta( helper ))
1297       ++nbAddedVols;
1298     else if ( nbNodes == 5 && _polygons.size() == 5 && addPyra ( helper ))
1299       ++nbAddedVols;
1300     else
1301     {
1302       ++nbAddedVols;
1303       helper.AddPolyhedralVolume( polyhedraNodes, quantities );
1304     }
1305     return nbAddedVols;
1306   }
1307
1308   //================================================================================
1309   /*!
1310    * \brief Return true if a polyhedron passes _sizeThreshold criterion
1311    */
1312   bool Hexahedron::checkPolyhedronSize() const
1313   {
1314     double volume = 0;
1315     for ( size_t iP = 0; iP < _polygons.size(); ++iP )
1316     {
1317       const _Face& polygon = _polygons[iP];
1318       gp_XYZ area (0,0,0);
1319       SMESH_TNodeXYZ p1 ( polygon._links[ 0 ].FirstNode()->Node() );
1320       for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
1321       {
1322         SMESH_TNodeXYZ p2 ( polygon._links[ iL ].LastNode()->Node() );
1323         area += p1 ^ p2;
1324         p1 = p2;
1325       }
1326       volume += p1 * area;
1327     }
1328     volume /= 6;
1329
1330     double initVolume = _sideLength[0] * _sideLength[1] * _sideLength[2];
1331
1332     return volume > initVolume / _sizeThreshold;
1333   }
1334   //================================================================================
1335   /*!
1336    * \brief Tries to create a hexahedron
1337    */
1338   bool Hexahedron::addHexa(SMESH_MesherHelper& helper)
1339   {
1340     if ( _polygons[0]._links.size() != 4 ||
1341          _polygons[1]._links.size() != 4 ||
1342          _polygons[2]._links.size() != 4 ||
1343          _polygons[3]._links.size() != 4 ||
1344          _polygons[4]._links.size() != 4 ||
1345          _polygons[5]._links.size() != 4   )
1346       return false;
1347     const SMDS_MeshNode* nodes[8];
1348     int nbN = 0;
1349     for ( int iL = 0; iL < 4; ++iL )
1350     {
1351       // a base node
1352       nodes[iL] = _polygons[0]._links[iL].FirstNode()->Node();
1353       ++nbN;
1354
1355       // find a top node above the base node
1356       _Link* link = _polygons[0]._links[iL]._link;
1357       ASSERT( link->_faces.size() > 1 );
1358       // a quadrangle sharing <link> with _polygons[0]
1359       _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[0] )];
1360       for ( int i = 0; i < 4; ++i )
1361         if ( quad->_links[i]._link == link )
1362         {
1363           // 1st node of a link opposite to <link> in <quad>
1364           nodes[iL+4] = quad->_links[(i+2)%4].FirstNode()->Node();
1365           ++nbN;
1366           break;
1367         }
1368     }
1369     if ( nbN == 8 )
1370       helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],
1371                         nodes[4],nodes[5],nodes[6],nodes[7] );
1372     return ( nbN == 8 );
1373   }
1374   //================================================================================
1375   /*!
1376    * \brief Tries to create a tetrahedron
1377    */
1378   bool Hexahedron::addTetra(SMESH_MesherHelper& helper)
1379   {
1380     const SMDS_MeshNode* nodes[4];
1381     nodes[0] = _polygons[0]._links[0].FirstNode()->Node();
1382     nodes[1] = _polygons[0]._links[1].FirstNode()->Node();
1383     nodes[2] = _polygons[0]._links[2].FirstNode()->Node();
1384
1385     _Link* link = _polygons[0]._links[0]._link;
1386     ASSERT( link->_faces.size() > 1 );
1387
1388     // a triangle sharing <link> with _polygons[0]
1389     _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[0] )];
1390     for ( int i = 0; i < 3; ++i )
1391       if ( tria->_links[i]._link == link )
1392       {
1393         nodes[3] = tria->_links[(i+1)%3].LastNode()->Node();
1394         helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3] );
1395         return true;
1396       }
1397
1398     return false;
1399   }
1400   //================================================================================
1401   /*!
1402    * \brief Tries to create a pentahedron
1403    */
1404   bool Hexahedron::addPenta(SMESH_MesherHelper& helper)
1405   {
1406     // find a base triangular face
1407     int iTri = -1;
1408     for ( int iF = 0; iF < 5 && iTri < 0; ++iF )
1409       if ( _polygons[ iF ]._links.size() == 3 )
1410         iTri = iF;
1411     if ( iTri < 0 ) return false;
1412
1413     // find nodes
1414     const SMDS_MeshNode* nodes[6];
1415     int nbN = 0;
1416     for ( int iL = 0; iL < 3; ++iL )
1417     {
1418       // a base node
1419       nodes[iL] = _polygons[ iTri ]._links[iL].FirstNode()->Node();
1420       ++nbN;
1421
1422       // find a top node above the base node
1423       _Link* link = _polygons[ iTri ]._links[iL]._link;
1424       ASSERT( link->_faces.size() > 1 );
1425       // a quadrangle sharing <link> with a base triangle
1426       _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[ iTri ] )];
1427       if ( quad->_links.size() != 4 ) return false;
1428       for ( int i = 0; i < 4; ++i )
1429         if ( quad->_links[i]._link == link )
1430         {
1431           // 1st node of a link opposite to <link> in <quad>
1432           nodes[iL+3] = quad->_links[(i+2)%4].FirstNode()->Node();
1433           ++nbN;
1434           break;
1435         }
1436     }
1437     if ( nbN == 6 )
1438       helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3], nodes[4],nodes[5] );
1439
1440     return ( nbN == 6 );
1441   }
1442   //================================================================================
1443   /*!
1444    * \brief Tries to create a pyramid
1445    */
1446   bool Hexahedron::addPyra(SMESH_MesherHelper& helper)
1447   {
1448     // find a base quadrangle
1449     int iQuad = -1;
1450     for ( int iF = 0; iF < 5 && iQuad < 0; ++iF )
1451       if ( _polygons[ iF ]._links.size() == 4 )
1452         iQuad = iF;
1453     if ( iQuad < 0 ) return false;
1454
1455     // find nodes
1456     const SMDS_MeshNode* nodes[5];
1457     nodes[0] = _polygons[iQuad]._links[0].FirstNode()->Node();
1458     nodes[1] = _polygons[iQuad]._links[1].FirstNode()->Node();
1459     nodes[2] = _polygons[iQuad]._links[2].FirstNode()->Node();
1460     nodes[3] = _polygons[iQuad]._links[3].FirstNode()->Node();
1461
1462     _Link* link = _polygons[iQuad]._links[0]._link;
1463     ASSERT( link->_faces.size() > 1 );
1464
1465     // a triangle sharing <link> with a base quadrangle
1466     _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[ iQuad ] )];
1467     if ( tria->_links.size() != 3 ) return false;
1468     for ( int i = 0; i < 3; ++i )
1469       if ( tria->_links[i]._link == link )
1470       {
1471         nodes[4] = tria->_links[(i+1)%3].LastNode()->Node();
1472         helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],nodes[4] );
1473         return true;
1474       }
1475
1476     return false;
1477   }
1478
1479 } // namespace
1480
1481 //=============================================================================
1482 /*!
1483  * \brief Generates 3D structured Cartesian mesh in the internal part of
1484  * solid shapes and polyhedral volumes near the shape boundary.
1485  *  \param theMesh - mesh to fill in
1486  *  \param theShape - a compound of all SOLIDs to mesh
1487  *  \retval bool - true in case of success
1488  */
1489 //=============================================================================
1490
1491 bool StdMeshers_Cartesian_3D::Compute(SMESH_Mesh &         theMesh,
1492                                       const TopoDS_Shape & theShape)
1493 {
1494   // The algorithm generates the mesh in following steps:
1495
1496   // 1) Intersection of grid lines with the geometry boundary.
1497   // This step allows to find out if a given node of the initial grid is
1498   // inside or outside the geometry.
1499
1500   // 2) For each cell of the grid, check how many of it's nodes are outside
1501   // of the geometry boundary. Depending on a result of this check
1502   // - skip a cell, if all it's nodes are outside
1503   // - skip a cell, if it is too small according to the size threshold
1504   // - add a hexahedron in the mesh, if all nodes are inside
1505   // - add a polyhedron in the mesh, if some nodes are inside and some outside
1506
1507   try {
1508     Grid grid;
1509
1510     TopTools_MapOfShape faceMap;
1511     for ( TopExp_Explorer fExp( theShape, TopAbs_FACE ); fExp.More(); fExp.Next() )
1512       if ( !faceMap.Add( fExp.Current() ))
1513         faceMap.Remove( fExp.Current() ); // remove a face shared by two solids
1514
1515     Bnd_Box shapeBox;
1516     vector<FaceGridIntersector> facesItersectors( faceMap.Extent() );
1517     TopTools_MapIteratorOfMapOfShape faceMppIt( faceMap );
1518     for ( int i = 0; faceMppIt.More(); faceMppIt.Next(), ++i )
1519     {
1520       facesItersectors[i]._face = TopoDS::Face( faceMppIt.Key() );
1521       facesItersectors[i]._grid = &grid;
1522       shapeBox.Add( facesItersectors[i].GetFaceBndBox() );
1523     }
1524
1525     vector<double> xCoords, yCoords, zCoords;
1526     _hyp->GetCoordinates( xCoords, yCoords, zCoords, shapeBox );
1527
1528     grid.SetCoordinates( xCoords, yCoords, zCoords, theShape );
1529
1530     // check if the grid encloses the shape
1531     if ( !_hyp->IsGridBySpacing(0) ||
1532          !_hyp->IsGridBySpacing(1) ||
1533          !_hyp->IsGridBySpacing(2) )
1534     {
1535       Bnd_Box gridBox;
1536       gridBox.Add( gp_Pnt( xCoords[0], yCoords[0], zCoords[0] ));
1537       gridBox.Add( gp_Pnt( xCoords.back(), yCoords.back(), zCoords.back() ));
1538       double x0,y0,z0, x1,y1,z1;
1539       shapeBox.Get(x0,y0,z0, x1,y1,z1);
1540       if ( gridBox.IsOut( gp_Pnt( x0,y0,z0 )) ||
1541            gridBox.IsOut( gp_Pnt( x1,y1,z1 )))
1542         for ( size_t i = 0; i < facesItersectors.size(); ++i )
1543           if ( !facesItersectors[i].IsInGrid( gridBox ))
1544             return error("The grid doesn't enclose the geometry");
1545     }
1546
1547     // Intersection of grid lines with the geometry boundary.
1548 #ifdef WITH_TBB
1549     tbb::parallel_for ( tbb::blocked_range<size_t>( 0, facesItersectors.size() ),
1550                         ParallelIntersector( facesItersectors ),
1551                         tbb::simple_partitioner());
1552 #else
1553     for ( size_t i = 0; i < facesItersectors.size(); ++i )
1554       facesItersectors[i].Intersect();
1555 #endif
1556
1557     // put interesection points onto the GridLine's; this is done after intersection
1558     // to avoid contention of facesItersectors for writing into the same GridLine
1559     // in case of parallel work of facesItersectors
1560     for ( size_t i = 0; i < facesItersectors.size(); ++i )
1561       facesItersectors[i].StoreIntersections();
1562
1563     SMESH_MesherHelper helper( theMesh );
1564     TopExp_Explorer solidExp (theShape, TopAbs_SOLID);
1565     helper.SetSubShape( solidExp.Current() );
1566     helper.SetElementsOnShape( true );
1567
1568     // create nodes on the geometry
1569     grid.ComputeNodes(helper);
1570
1571     Hexahedron hex( _hyp->GetSizeThreshold(), &grid );
1572     int nbAdded = 0;
1573     for ( size_t k = 1; k < zCoords.size(); ++k )
1574       for ( size_t j = 1; j < yCoords.size(); ++j )
1575         for ( size_t i = 1; i < xCoords.size(); ++i )
1576         {
1577           hex.Init( i-1, j-1, k-1 );
1578           nbAdded += hex.MakeElements( helper );
1579         }
1580
1581     SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1582     if ( nbAdded > 0 )
1583     {
1584       // make all SOLIDS computed
1585       if ( SMESHDS_SubMesh* sm1 = meshDS->MeshElements( solidExp.Current()) )
1586       {
1587         SMDS_ElemIteratorPtr volIt = sm1->GetElements();
1588         for ( ; solidExp.More() && volIt->more(); solidExp.Next() )
1589         {
1590           const SMDS_MeshElement* vol = volIt->next();
1591           sm1->RemoveElement( vol, /*isElemDeleted=*/false );
1592           meshDS->SetMeshElementOnShape( vol, solidExp.Current() );
1593         }
1594       }
1595       // make other sub-shapes computed
1596       setSubmeshesComputed( theMesh, theShape );
1597       return true;
1598     }
1599
1600     // remove free nodes
1601     if ( SMESHDS_SubMesh * smDS = meshDS->MeshElements( helper.GetSubShapeID() ))
1602     {
1603       // grid nodes
1604       for ( size_t i = 0; i < grid._nodes.size(); ++i )
1605         if ( grid._nodes[i] && grid._nodes[i]->NbInverseElements() == 0 )
1606           meshDS->RemoveFreeNode( grid._nodes[i], smDS, /*fromGroups=*/false );
1607
1608       // intersection nodes
1609       for ( int iDir = 0; iDir < 3; ++iDir )
1610       {
1611         vector< GridLine >& lines = grid._lines[ iDir ];
1612         for ( size_t i = 0; i < lines.size(); ++i )
1613         {
1614           multiset< IntersectionPoint >::iterator ip = lines[i]._intPoints.begin();
1615           for ( ; ip != lines[i]._intPoints.end(); ++ip )
1616             if ( ip->_node && ip->_node->NbInverseElements() == 0 )
1617               meshDS->RemoveFreeNode( grid._nodes[i], smDS, /*fromGroups=*/false );
1618         }
1619       }
1620     }
1621
1622     // TODO: evalute time
1623
1624   }
1625   // SMESH_ComputeError is not caught at SMESH_submesh level for an unknown reason
1626   catch ( SMESH_ComputeError& e)
1627   {
1628     return error( SMESH_ComputeErrorPtr( new SMESH_ComputeError( e )));
1629   }
1630   return false;
1631 }
1632
1633 //=============================================================================
1634 /*!
1635  *  Evaluate
1636  */
1637 //=============================================================================
1638
1639 bool StdMeshers_Cartesian_3D::Evaluate(SMESH_Mesh &         theMesh,
1640                                        const TopoDS_Shape & theShape,
1641                                        MapShapeNbElems&     theResMap)
1642 {
1643   // TODO
1644 //   std::vector<int> aResVec(SMDSEntity_Last);
1645 //   for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
1646 //   if(IsQuadratic) {
1647 //     aResVec[SMDSEntity_Quad_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
1648 //     int nb1d_face0_int = ( nb2d_face0*4 - nb1d ) / 2;
1649 //     aResVec[SMDSEntity_Node] = nb0d_face0 * ( 2*nb2d/nb1d - 1 ) - nb1d_face0_int * nb2d/nb1d;
1650 //   }
1651 //   else {
1652 //     aResVec[SMDSEntity_Node] = nb0d_face0 * ( nb2d/nb1d - 1 );
1653 //     aResVec[SMDSEntity_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
1654 //   }
1655 //   SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
1656 //   aResMap.insert(std::make_pair(sm,aResVec));
1657
1658   return true;
1659 }
1660
1661 //=============================================================================
1662 namespace
1663 {
1664   /*!
1665    * \brief Event listener setting/unsetting _alwaysComputed flag to
1666    *        submeshes of inferior levels to avoid their computing
1667    */
1668   struct _EventListener : public SMESH_subMeshEventListener
1669   {
1670     string _algoName;
1671
1672     _EventListener(const string& algoName):
1673       SMESH_subMeshEventListener(/*isDeletable=*/true), _algoName(algoName) {}
1674
1675     // --------------------------------------------------------------------------------
1676     // setting/unsetting _alwaysComputed flag to submeshes of inferior levels
1677     //
1678     static void setAlwaysComputed( const bool     isComputed,
1679                                    SMESH_subMesh* subMeshOfSolid)
1680     {
1681       SMESH_subMeshIteratorPtr smIt =
1682         subMeshOfSolid->getDependsOnIterator(/*includeSelf=*/false, /*complexShapeFirst=*/false);
1683       while ( smIt->more() )
1684       {
1685         SMESH_subMesh* sm = smIt->next();
1686         sm->SetIsAlwaysComputed( isComputed );
1687       }
1688     }
1689
1690     // --------------------------------------------------------------------------------
1691     // unsetting _alwaysComputed flag if "Cartesian_3D" was removed
1692     //
1693     virtual void ProcessEvent(const int          event,
1694                               const int          eventType,
1695                               SMESH_subMesh*     subMeshOfSolid,
1696                               SMESH_subMeshEventListenerData* data,
1697                               const SMESH_Hypothesis*         hyp = 0)
1698     {
1699       if ( eventType == SMESH_subMesh::COMPUTE_EVENT )
1700       {
1701         setAlwaysComputed( subMeshOfSolid->GetComputeState() == SMESH_subMesh::COMPUTE_OK,
1702                            subMeshOfSolid );
1703       }
1704       else
1705       {
1706         SMESH_Algo* algo3D = subMeshOfSolid->GetAlgo();
1707         if ( !algo3D || _algoName != algo3D->GetName() )
1708           setAlwaysComputed( false, subMeshOfSolid );
1709       }
1710     }
1711
1712     // --------------------------------------------------------------------------------
1713     // set the event listener
1714     //
1715     static void SetOn( SMESH_subMesh* subMeshOfSolid, const string& algoName )
1716     {
1717       subMeshOfSolid->SetEventListener( new _EventListener( algoName ),
1718                                         /*data=*/0,
1719                                         subMeshOfSolid );
1720     }
1721
1722   }; // struct _EventListener
1723
1724 } // namespace
1725
1726 //================================================================================
1727 /*!
1728  * \brief Sets event listener to submeshes if necessary
1729  *  \param subMesh - submesh where algo is set
1730  * This method is called when a submesh gets HYP_OK algo_state.
1731  * After being set, event listener is notified on each event of a submesh.
1732  */
1733 //================================================================================
1734
1735 void StdMeshers_Cartesian_3D::SetEventListener(SMESH_subMesh* subMesh)
1736 {
1737   _EventListener::SetOn( subMesh, GetName() );
1738 }
1739
1740 //================================================================================
1741 /*!
1742  * \brief Set _alwaysComputed flag to submeshes of inferior levels to avoid their computing
1743  */
1744 //================================================================================
1745
1746 void StdMeshers_Cartesian_3D::setSubmeshesComputed(SMESH_Mesh&         theMesh,
1747                                                    const TopoDS_Shape& theShape)
1748 {
1749   for ( TopExp_Explorer soExp( theShape, TopAbs_SOLID ); soExp.More(); soExp.Next() )
1750     _EventListener::setAlwaysComputed( true, theMesh.GetSubMesh( soExp.Current() ));
1751 }
1752