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