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