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