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