Salome HOME
Correct python dump for the SHAPERSTUDY part of the script
[modules/smesh.git] / src / StdMeshers / StdMeshers_Cartesian_3D.cxx
1 // Copyright (C) 2007-2019  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, or (at your option) any later version.
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 "SMESHDS_Mesh.hxx"
29 #include "SMESH_Block.hxx"
30 #include "SMESH_Comment.hxx"
31 #include "SMESH_Mesh.hxx"
32 #include "SMESH_MesherHelper.hxx"
33 #include "SMESH_subMesh.hxx"
34 #include "SMESH_subMeshEventListener.hxx"
35 #include "StdMeshers_CartesianParameters3D.hxx"
36
37 #include <utilities.h>
38 #include <Utils_ExceptHandlers.hxx>
39
40 #include <GEOMUtils.hxx>
41
42 #include <BRepAdaptor_Curve.hxx>
43 #include <BRepAdaptor_Surface.hxx>
44 #include <BRepBndLib.hxx>
45 #include <BRepBuilderAPI_Copy.hxx>
46 #include <BRepBuilderAPI_MakeFace.hxx>
47 #include <BRepTools.hxx>
48 #include <BRepTopAdaptor_FClass2d.hxx>
49 #include <BRep_Builder.hxx>
50 #include <BRep_Tool.hxx>
51 #include <Bnd_B3d.hxx>
52 #include <Bnd_Box.hxx>
53 #include <ElSLib.hxx>
54 #include <GCPnts_UniformDeflection.hxx>
55 #include <Geom2d_BSplineCurve.hxx>
56 #include <Geom2d_BezierCurve.hxx>
57 #include <Geom2d_TrimmedCurve.hxx>
58 #include <GeomAPI_ProjectPointOnSurf.hxx>
59 #include <GeomLib.hxx>
60 #include <Geom_BSplineCurve.hxx>
61 #include <Geom_BSplineSurface.hxx>
62 #include <Geom_BezierCurve.hxx>
63 #include <Geom_BezierSurface.hxx>
64 #include <Geom_RectangularTrimmedSurface.hxx>
65 #include <Geom_TrimmedCurve.hxx>
66 #include <IntAna_IntConicQuad.hxx>
67 #include <IntAna_IntLinTorus.hxx>
68 #include <IntAna_Quadric.hxx>
69 #include <IntCurveSurface_TransitionOnCurve.hxx>
70 #include <IntCurvesFace_Intersector.hxx>
71 #include <Poly_Triangulation.hxx>
72 #include <Precision.hxx>
73 #include <TopExp.hxx>
74 #include <TopExp_Explorer.hxx>
75 #include <TopLoc_Location.hxx>
76 #include <TopTools_IndexedMapOfShape.hxx>
77 #include <TopTools_MapOfShape.hxx>
78 #include <TopoDS.hxx>
79 #include <TopoDS_Compound.hxx>
80 #include <TopoDS_Face.hxx>
81 #include <TopoDS_TShape.hxx>
82 #include <gp_Cone.hxx>
83 #include <gp_Cylinder.hxx>
84 #include <gp_Lin.hxx>
85 #include <gp_Pln.hxx>
86 #include <gp_Pnt2d.hxx>
87 #include <gp_Sphere.hxx>
88 #include <gp_Torus.hxx>
89
90 #include <limits>
91
92 //#undef WITH_TBB
93 #ifdef WITH_TBB
94
95 #ifdef WIN32
96 // See https://docs.microsoft.com/en-gb/cpp/porting/modifying-winver-and-win32-winnt?view=vs-2019
97 // Windows 10 = 0x0A00  
98 #define WINVER 0x0A00
99 #define _WIN32_WINNT 0x0A00
100
101 #endif
102
103 #include <tbb/parallel_for.h>
104 //#include <tbb/enumerable_thread_specific.h>
105 #endif
106
107 using namespace std;
108
109 #ifdef _DEBUG_
110 //#define _MY_DEBUG_
111 #endif
112
113 //=============================================================================
114 /*!
115  * Constructor
116  */
117 //=============================================================================
118
119 StdMeshers_Cartesian_3D::StdMeshers_Cartesian_3D(int hypId, SMESH_Gen * gen)
120   :SMESH_3D_Algo(hypId, gen)
121 {
122   _name = "Cartesian_3D";
123   _shapeType = (1 << TopAbs_SOLID);       // 1 bit /shape type
124   _compatibleHypothesis.push_back("CartesianParameters3D");
125
126   _onlyUnaryInput = false;          // to mesh all SOLIDs at once
127   _requireDiscreteBoundary = false; // 2D mesh not needed
128   _supportSubmeshes = false;        // do not use any existing mesh
129 }
130
131 //=============================================================================
132 /*!
133  * Check presence of a hypothesis
134  */
135 //=============================================================================
136
137 bool StdMeshers_Cartesian_3D::CheckHypothesis (SMESH_Mesh&          aMesh,
138                                                const TopoDS_Shape&  aShape,
139                                                Hypothesis_Status&   aStatus)
140 {
141   aStatus = SMESH_Hypothesis::HYP_MISSING;
142
143   const list<const SMESHDS_Hypothesis*>& hyps = GetUsedHypothesis(aMesh, aShape);
144   list <const SMESHDS_Hypothesis* >::const_iterator h = hyps.begin();
145   if ( h == hyps.end())
146   {
147     return false;
148   }
149
150   for ( ; h != hyps.end(); ++h )
151   {
152     if (( _hyp = dynamic_cast<const StdMeshers_CartesianParameters3D*>( *h )))
153     {
154       aStatus = _hyp->IsDefined() ? HYP_OK : HYP_BAD_PARAMETER;
155       break;
156     }
157   }
158
159   return aStatus == HYP_OK;
160 }
161
162 namespace
163 {
164   typedef int TGeomID;
165
166   //=============================================================================
167   // Definitions of internal utils
168   // --------------------------------------------------------------------------
169   enum Transition {
170     Trans_TANGENT = IntCurveSurface_Tangent,
171     Trans_IN      = IntCurveSurface_In,
172     Trans_OUT     = IntCurveSurface_Out,
173     Trans_APEX
174   };
175   // --------------------------------------------------------------------------
176   /*!
177    * \brief Common data of any intersection between a Grid and a shape
178    */
179   struct B_IntersectPoint
180   {
181     mutable const SMDS_MeshNode* _node;
182     mutable vector< TGeomID >    _faceIDs;
183
184     B_IntersectPoint(): _node(NULL) {}
185     void Add( const vector< TGeomID >& fIDs, const SMDS_MeshNode* n=0 ) const;
186     int HasCommonFace( const B_IntersectPoint * other, int avoidFace=-1 ) const;
187     bool IsOnFace( int faceID ) const;
188     virtual ~B_IntersectPoint() {}
189   };
190   // --------------------------------------------------------------------------
191   /*!
192    * \brief Data of intersection between a GridLine and a TopoDS_Face
193    */
194   struct F_IntersectPoint : public B_IntersectPoint
195   {
196     double             _paramOnLine;
197     mutable Transition _transition;
198     mutable size_t     _indexOnLine;
199
200     bool operator< ( const F_IntersectPoint& o ) const { return _paramOnLine < o._paramOnLine; }
201   };
202   // --------------------------------------------------------------------------
203   /*!
204    * \brief Data of intersection between GridPlanes and a TopoDS_EDGE
205    */
206   struct E_IntersectPoint : public B_IntersectPoint
207   {
208     gp_Pnt  _point;
209     double  _uvw[3];
210     TGeomID _shapeID;
211   };
212   // --------------------------------------------------------------------------
213   /*!
214    * \brief A line of the grid and its intersections with 2D geometry
215    */
216   struct GridLine
217   {
218     gp_Lin _line;
219     double _length; // line length
220     multiset< F_IntersectPoint > _intPoints;
221
222     void RemoveExcessIntPoints( const double tol );
223     bool GetIsOutBefore( multiset< F_IntersectPoint >::iterator ip, bool prevIsOut );
224   };
225   // --------------------------------------------------------------------------
226   /*!
227    * \brief Planes of the grid used to find intersections of an EDGE with a hexahedron
228    */
229   struct GridPlanes
230   {
231     gp_XYZ           _zNorm;
232     vector< gp_XYZ > _origins; // origin points of all planes in one direction
233     vector< double > _zProjs;  // projections of origins to _zNorm
234   };
235   // --------------------------------------------------------------------------
236   /*!
237    * \brief Iterator on the parallel grid lines of one direction
238    */
239   struct LineIndexer
240   {
241     size_t _size  [3];
242     size_t _curInd[3];
243     size_t _iVar1, _iVar2, _iConst;
244     string _name1, _name2, _nameConst;
245     LineIndexer() {}
246     LineIndexer( size_t sz1, size_t sz2, size_t sz3,
247                  size_t iv1, size_t iv2, size_t iConst,
248                  const string& nv1, const string& nv2, const string& nConst )
249     {
250       _size[0] = sz1; _size[1] = sz2; _size[2] = sz3;
251       _curInd[0] = _curInd[1] = _curInd[2] = 0;
252       _iVar1 = iv1; _iVar2 = iv2; _iConst = iConst; 
253       _name1 = nv1; _name2 = nv2; _nameConst = nConst;
254     }
255
256     size_t I() const { return _curInd[0]; }
257     size_t J() const { return _curInd[1]; }
258     size_t K() const { return _curInd[2]; }
259     void SetIJK( size_t i, size_t j, size_t k )
260     {
261       _curInd[0] = i; _curInd[1] = j; _curInd[2] = k;
262     }
263     void operator++()
264     {
265       if ( ++_curInd[_iVar1] == _size[_iVar1] )
266         _curInd[_iVar1] = 0, ++_curInd[_iVar2];
267     }
268     bool More() const { return _curInd[_iVar2] < _size[_iVar2]; }
269     size_t LineIndex   () const { return _curInd[_iVar1] + _curInd[_iVar2]* _size[_iVar1]; }
270     size_t LineIndex10 () const { return (_curInd[_iVar1] + 1 ) + _curInd[_iVar2]* _size[_iVar1]; }
271     size_t LineIndex01 () const { return _curInd[_iVar1] + (_curInd[_iVar2] + 1 )* _size[_iVar1]; }
272     size_t LineIndex11 () const { return (_curInd[_iVar1] + 1 ) + (_curInd[_iVar2] + 1 )* _size[_iVar1]; }
273     void SetIndexOnLine (size_t i)  { _curInd[ _iConst ] = i; }
274     size_t NbLines() const { return _size[_iVar1] * _size[_iVar2]; }
275   };
276   // --------------------------------------------------------------------------
277   /*!
278    * \brief Container of GridLine's
279    */
280   struct Grid
281   {
282     vector< double >   _coords[3]; // coordinates of grid nodes
283     gp_XYZ             _axes  [3]; // axis directions
284     vector< GridLine > _lines [3]; //    in 3 directions
285     double             _tol, _minCellSize;
286     gp_XYZ             _origin;
287     gp_Mat             _invB; // inverted basis of _axes
288
289     vector< const SMDS_MeshNode* >    _nodes; // mesh nodes at grid nodes
290     vector< const F_IntersectPoint* > _gridIntP; // grid node intersection with geometry
291
292     list< E_IntersectPoint >          _edgeIntP; // intersections with EDGEs
293     TopTools_IndexedMapOfShape        _shapes;
294
295     SMESH_MesherHelper*               _helper;
296
297     size_t CellIndex( size_t i, size_t j, size_t k ) const
298     {
299       return i + j*(_coords[0].size()-1) + k*(_coords[0].size()-1)*(_coords[1].size()-1);
300     }
301     size_t NodeIndex( size_t i, size_t j, size_t k ) const
302     {
303       return i + j*_coords[0].size() + k*_coords[0].size()*_coords[1].size();
304     }
305     size_t NodeIndexDX() const { return 1; }
306     size_t NodeIndexDY() const { return _coords[0].size(); }
307     size_t NodeIndexDZ() const { return _coords[0].size() * _coords[1].size(); }
308
309     LineIndexer GetLineIndexer(size_t iDir) const;
310
311     void SetCoordinates(const vector<double>& xCoords,
312                         const vector<double>& yCoords,
313                         const vector<double>& zCoords,
314                         const double*         axesDirs,
315                         const Bnd_Box&        bndBox );
316     void ComputeUVW(const gp_XYZ& p, double uvw[3]);
317     void ComputeNodes(SMESH_MesherHelper& helper);
318   };
319   // --------------------------------------------------------------------------
320   /*!
321    * \brief Intersector of TopoDS_Face with all GridLine's
322    */
323   struct FaceGridIntersector
324   {
325     TopoDS_Face _face;
326     TGeomID     _faceID;
327     Grid*       _grid;
328     Bnd_Box     _bndBox;
329     IntCurvesFace_Intersector* _surfaceInt;
330     vector< std::pair< GridLine*, F_IntersectPoint > > _intersections;
331
332     FaceGridIntersector(): _grid(0), _surfaceInt(0) {}
333     void Intersect();
334
335     void StoreIntersections()
336     {
337       for ( size_t i = 0; i < _intersections.size(); ++i )
338       {
339         multiset< F_IntersectPoint >::iterator ip = 
340           _intersections[i].first->_intPoints.insert( _intersections[i].second );
341         ip->_faceIDs.reserve( 1 );
342         ip->_faceIDs.push_back( _faceID );
343       }
344     }
345     const Bnd_Box& GetFaceBndBox()
346     {
347       GetCurveFaceIntersector();
348       return _bndBox;
349     }
350     IntCurvesFace_Intersector* GetCurveFaceIntersector()
351     {
352       if ( !_surfaceInt )
353       {
354         _surfaceInt = new IntCurvesFace_Intersector( _face, Precision::PConfusion() );
355         _bndBox     = _surfaceInt->Bounding();
356         if ( _bndBox.IsVoid() )
357           BRepBndLib::Add (_face, _bndBox);
358       }
359       return _surfaceInt;
360     }
361     bool IsThreadSafe(set< const Standard_Transient* >& noSafeTShapes) const;
362   };
363   // --------------------------------------------------------------------------
364   /*!
365    * \brief Intersector of a surface with a GridLine
366    */
367   struct FaceLineIntersector
368   {
369     double      _tol;
370     double      _u, _v, _w; // params on the face and the line
371     Transition  _transition; // transition of at intersection (see IntCurveSurface.cdl)
372     Transition  _transIn, _transOut; // IN and OUT transitions depending of face orientation
373
374     gp_Pln      _plane;
375     gp_Cylinder _cylinder;
376     gp_Cone     _cone;
377     gp_Sphere   _sphere;
378     gp_Torus    _torus;
379     IntCurvesFace_Intersector* _surfaceInt;
380
381     vector< F_IntersectPoint > _intPoints;
382
383     void IntersectWithPlane   (const GridLine& gridLine);
384     void IntersectWithCylinder(const GridLine& gridLine);
385     void IntersectWithCone    (const GridLine& gridLine);
386     void IntersectWithSphere  (const GridLine& gridLine);
387     void IntersectWithTorus   (const GridLine& gridLine);
388     void IntersectWithSurface (const GridLine& gridLine);
389
390     bool UVIsOnFace() const;
391     void addIntPoint(const bool toClassify=true);
392     bool isParamOnLineOK( const double linLength )
393     {
394       return -_tol < _w && _w < linLength + _tol;
395     }
396     FaceLineIntersector():_surfaceInt(0) {}
397     ~FaceLineIntersector() { if (_surfaceInt ) delete _surfaceInt; _surfaceInt = 0; }
398   };
399   // --------------------------------------------------------------------------
400   /*!
401    * \brief Class representing topology of the hexahedron and creating a mesh
402    *        volume basing on analysis of hexahedron intersection with geometry
403    */
404   class Hexahedron
405   {
406     // --------------------------------------------------------------------------------
407     struct _Face;
408     struct _Link;
409     // --------------------------------------------------------------------------------
410     struct _Node //!< node either at a hexahedron corner or at intersection
411     {
412       const SMDS_MeshNode*    _node; // mesh node at hexahedron corner
413       const B_IntersectPoint* _intPoint;
414       const _Face*            _usedInFace;
415
416       _Node(const SMDS_MeshNode* n=0, const B_IntersectPoint* ip=0)
417         :_node(n), _intPoint(ip), _usedInFace(0) {} 
418       const SMDS_MeshNode*    Node() const
419       { return ( _intPoint && _intPoint->_node ) ? _intPoint->_node : _node; }
420       const E_IntersectPoint* EdgeIntPnt() const
421       { return static_cast< const E_IntersectPoint* >( _intPoint ); }
422       bool IsUsedInFace( const _Face* polygon = 0 )
423       {
424         return polygon ? ( _usedInFace == polygon ) : bool( _usedInFace );
425       }
426       void Add( const E_IntersectPoint* ip )
427       {
428         if ( !_intPoint ) {
429           _intPoint = ip;
430         }
431         else if ( !_intPoint->_node ) {
432           ip->Add( _intPoint->_faceIDs );
433           _intPoint = ip;
434         }
435         else  {
436           _intPoint->Add( ip->_faceIDs );
437         }
438       }
439       TGeomID IsLinked( const B_IntersectPoint* other,
440                         TGeomID                 avoidFace=-1 ) const // returns id of a common face
441       {
442         return _intPoint ? _intPoint->HasCommonFace( other, avoidFace ) : 0;
443       }
444       bool IsOnFace( TGeomID faceID ) const // returns true if faceID is found
445       {
446         return _intPoint ? _intPoint->IsOnFace( faceID ) : false;
447       }
448       gp_Pnt Point() const
449       {
450         if ( const SMDS_MeshNode* n = Node() )
451           return SMESH_TNodeXYZ( n );
452         if ( const E_IntersectPoint* eip =
453              dynamic_cast< const E_IntersectPoint* >( _intPoint ))
454           return eip->_point;
455         return gp_Pnt( 1e100, 0, 0 );
456       }
457       TGeomID ShapeID() const
458       {
459         if ( const E_IntersectPoint* eip = dynamic_cast< const E_IntersectPoint* >( _intPoint ))
460           return eip->_shapeID;
461         return 0;
462       }
463     };
464     // --------------------------------------------------------------------------------
465     struct _Link // link connecting two _Node's
466     {
467       _Node* _nodes[2];
468       _Face* _faces[2]; // polygons sharing a link
469       vector< const F_IntersectPoint* > _fIntPoints; // GridLine intersections with FACEs
470       vector< _Node* >                  _fIntNodes;   // _Node's at _fIntPoints
471       vector< _Link >                   _splits;
472       _Link() { _faces[0] = 0; }
473     };
474     // --------------------------------------------------------------------------------
475     struct _OrientedLink
476     {
477       _Link* _link;
478       bool   _reverse;
479       _OrientedLink( _Link* link=0, bool reverse=false ): _link(link), _reverse(reverse) {}
480       void Reverse() { _reverse = !_reverse; }
481       int NbResultLinks() const { return _link->_splits.size(); }
482       _OrientedLink ResultLink(int i) const
483       {
484         return _OrientedLink(&_link->_splits[_reverse ? NbResultLinks()-i-1 : i],_reverse);
485       }
486       _Node* FirstNode() const { return _link->_nodes[ _reverse ]; }
487       _Node* LastNode()  const { return _link->_nodes[ !_reverse ]; }
488       operator bool() const { return _link; }
489       vector< TGeomID > GetNotUsedFace(const set<TGeomID>& usedIDs ) const // returns supporting FACEs
490       {
491         vector< TGeomID > faces;
492         const B_IntersectPoint *ip0, *ip1;
493         if (( ip0 = _link->_nodes[0]->_intPoint ) &&
494             ( ip1 = _link->_nodes[1]->_intPoint ))
495         {
496           for ( size_t i = 0; i < ip0->_faceIDs.size(); ++i )
497             if ( ip1->IsOnFace ( ip0->_faceIDs[i] ) &&
498                  !usedIDs.count( ip0->_faceIDs[i] ) )
499               faces.push_back( ip0->_faceIDs[i] );
500         }
501         return faces;
502       }
503       bool HasEdgeNodes() const
504       {
505         return ( dynamic_cast< const E_IntersectPoint* >( _link->_nodes[0]->_intPoint ) ||
506                  dynamic_cast< const E_IntersectPoint* >( _link->_nodes[1]->_intPoint ));
507       }
508       int NbFaces() const
509       {
510         return !_link->_faces[0] ? 0 : 1 + bool( _link->_faces[1] );
511       }
512       void AddFace( _Face* f )
513       {
514         if ( _link->_faces[0] )
515         {
516           _link->_faces[1] = f;
517         }
518         else
519         {
520           _link->_faces[0] = f;
521           _link->_faces[1] = 0;
522         }
523       }
524       void RemoveFace( _Face* f )
525       {
526         if ( !_link->_faces[0] ) return;
527
528         if ( _link->_faces[1] == f )
529         {
530           _link->_faces[1] = 0;
531         }
532         else if ( _link->_faces[0] == f )
533         {
534           _link->_faces[0] = 0;
535           if ( _link->_faces[1] )
536           {
537             _link->_faces[0] = _link->_faces[1];
538             _link->_faces[1] = 0;
539           }
540         }
541       }
542     };
543     // --------------------------------------------------------------------------------
544     struct _Face
545     {
546       vector< _OrientedLink > _links;       // links on GridLine's
547       vector< _Link >         _polyLinks;   // links added to close a polygonal face
548       vector< _Node* >        _eIntNodes;   // nodes at intersection with EDGEs
549       bool IsPolyLink( const _OrientedLink& ol )
550       {
551         return _polyLinks.empty() ? false :
552           ( &_polyLinks[0] <= ol._link &&  ol._link <= &_polyLinks.back() );
553       }
554       void AddPolyLink(_Node* n0, _Node* n1, _Face* faceToFindEqual=0)
555       {
556         if ( faceToFindEqual && faceToFindEqual != this ) {
557           for ( size_t iL = 0; iL < faceToFindEqual->_polyLinks.size(); ++iL )
558             if ( faceToFindEqual->_polyLinks[iL]._nodes[0] == n1 &&
559                  faceToFindEqual->_polyLinks[iL]._nodes[1] == n0 )
560             {
561               _links.push_back
562                 ( _OrientedLink( & faceToFindEqual->_polyLinks[iL], /*reverse=*/true ));
563               return;
564             }
565         }
566         _Link l;
567         l._nodes[0] = n0;
568         l._nodes[1] = n1;
569         _polyLinks.push_back( l );
570         _links.push_back( _OrientedLink( &_polyLinks.back() ));
571       }
572     };
573     // --------------------------------------------------------------------------------
574     struct _volumeDef // holder of nodes of a volume mesh element
575     {
576       vector< _Node* > _nodes;
577       vector< int >    _quantities;
578       typedef boost::shared_ptr<_volumeDef> Ptr;
579       void set( const vector< _Node* >& nodes,
580                 const vector< int >&    quant = vector< int >() )
581       { _nodes = nodes; _quantities = quant; }
582       void set( _Node** nodes, int nb )
583       { _nodes.assign( nodes, nodes + nb ); }
584     };
585
586     // topology of a hexahedron
587     int   _nodeShift[8];
588     _Node _hexNodes [8];
589     _Link _hexLinks [12];
590     _Face _hexQuads [6];
591
592     // faces resulted from hexahedron intersection
593     vector< _Face > _polygons;
594
595     // intresections with EDGEs
596     vector< const E_IntersectPoint* > _eIntPoints;
597
598     // additional nodes created at intersection points
599     vector< _Node > _intNodes;
600
601     // nodes inside the hexahedron (at VERTEXes)
602     vector< _Node* > _vIntNodes;
603
604     // computed volume elements
605     //vector< _volumeDef::Ptr > _volumeDefs;
606     _volumeDef _volumeDefs;
607
608     Grid*       _grid;
609     double      _sizeThreshold, _sideLength[3];
610     int         _nbCornerNodes, _nbFaceIntNodes, _nbBndNodes;
611     int         _origNodeInd; // index of _hexNodes[0] node within the _grid
612     size_t      _i,_j,_k;
613
614   public:
615     Hexahedron(const double sizeThreshold, Grid* grid);
616     int MakeElements(SMESH_MesherHelper&                      helper,
617                      const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap);
618     void ComputeElements();
619     void Init() { init( _i, _j, _k ); }
620
621   private:
622     Hexahedron(const Hexahedron& other );
623     void init( size_t i, size_t j, size_t k );
624     void init( size_t i );
625     void addEdges(SMESH_MesherHelper&                      helper,
626                   vector< Hexahedron* >&                   intersectedHex,
627                   const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap);
628     gp_Pnt findIntPoint( double u1, double proj1, double u2, double proj2,
629                          double proj, BRepAdaptor_Curve& curve,
630                          const gp_XYZ& axis, const gp_XYZ& origin );
631     int  getEntity( const E_IntersectPoint* ip, int* facets, int& sub );
632     bool addIntersection( const E_IntersectPoint& ip,
633                           vector< Hexahedron* >&  hexes,
634                           int ijk[], int dIJK[] );
635     bool findChain( _Node* n1, _Node* n2, _Face& quad, vector<_Node*>& chainNodes );
636     bool closePolygon( _Face* polygon, vector<_Node*>& chainNodes ) const;
637     bool findChainOnEdge( const vector< _OrientedLink >& splits,
638                           const _OrientedLink&           prevSplit,
639                           const _OrientedLink&           avoidSplit,
640                           size_t &                       iS,
641                           _Face&                         quad,
642                           vector<_Node*>&                chn);
643     int  addElements(SMESH_MesherHelper& helper);
644     bool isOutPoint( _Link& link, int iP, SMESH_MesherHelper& helper ) const;
645     void sortVertexNodes(vector<_Node*>& nodes, _Node* curNode, TGeomID face);
646     bool isInHole() const;
647     bool checkPolyhedronSize() const;
648     bool addHexa ();
649     bool addTetra();
650     bool addPenta();
651     bool addPyra ();
652     bool debugDumpLink( _Link* link );
653     _Node* findEqualNode( vector< _Node* >&       nodes,
654                           const E_IntersectPoint* ip,
655                           const double            tol2 )
656     {
657       for ( size_t i = 0; i < nodes.size(); ++i )
658         if ( nodes[i]->EdgeIntPnt() == ip ||
659              nodes[i]->Point().SquareDistance( ip->_point ) <= tol2 )
660           return nodes[i];
661       return 0;
662     }
663     bool isImplementEdges() const { return !_grid->_edgeIntP.empty(); }
664     bool isOutParam(const double uvw[3]) const;
665   };
666
667 #ifdef WITH_TBB
668   // --------------------------------------------------------------------------
669   /*!
670    * \brief Hexahedron computing volumes in one thread
671    */
672   struct ParallelHexahedron
673   {
674     vector< Hexahedron* >& _hexVec;
675     ParallelHexahedron( vector< Hexahedron* >& hv ): _hexVec(hv) {}
676     void operator() ( const tbb::blocked_range<size_t>& r ) const
677     {
678       for ( size_t i = r.begin(); i != r.end(); ++i )
679         if ( Hexahedron* hex = _hexVec[ i ] )
680           hex->ComputeElements();
681     }
682   };
683   // --------------------------------------------------------------------------
684   /*!
685    * \brief Structure intersecting certain nb of faces with GridLine's in one thread
686    */
687   struct ParallelIntersector
688   {
689     vector< FaceGridIntersector >& _faceVec;
690     ParallelIntersector( vector< FaceGridIntersector >& faceVec): _faceVec(faceVec){}
691     void operator() ( const tbb::blocked_range<size_t>& r ) const
692     {
693       for ( size_t i = r.begin(); i != r.end(); ++i )
694         _faceVec[i].Intersect();
695     }
696   };
697 #endif
698
699   //=============================================================================
700   // Implementation of internal utils
701   //=============================================================================
702   /*!
703    * \brief adjust \a i to have \a val between values[i] and values[i+1]
704    */
705   inline void locateValue( int & i, double val, const vector<double>& values,
706                            int& di, double tol )
707   {
708     //val += values[0]; // input \a val is measured from 0.
709     if ( i > (int) values.size()-2 )
710       i = values.size()-2;
711     else
712       while ( i+2 < (int) values.size() && val > values[ i+1 ])
713         ++i;
714     while ( i > 0 && val < values[ i ])
715       --i;
716
717     if ( i > 0 && val - values[ i ] < tol )
718       di = -1;
719     else if ( i+2 < (int) values.size() && values[ i+1 ] - val < tol )
720       di = 1;
721     else
722       di = 0;
723   }
724   //=============================================================================
725   /*
726    * Remove coincident intersection points
727    */
728   void GridLine::RemoveExcessIntPoints( const double tol )
729   {
730     if ( _intPoints.size() < 2 ) return;
731
732     set< Transition > tranSet;
733     multiset< F_IntersectPoint >::iterator ip1, ip2 = _intPoints.begin();
734     while ( ip2 != _intPoints.end() )
735     {
736       tranSet.clear();
737       ip1 = ip2++;
738       while ( ip2 != _intPoints.end() && ip2->_paramOnLine - ip1->_paramOnLine <= tol )
739       {
740         tranSet.insert( ip1->_transition );
741         tranSet.insert( ip2->_transition );
742         ip2->Add( ip1->_faceIDs );
743         _intPoints.erase( ip1 );
744         ip1 = ip2++;
745       }
746       if ( tranSet.size() > 1 ) // points with different transition coincide
747       {
748         bool isIN  = tranSet.count( Trans_IN );
749         bool isOUT = tranSet.count( Trans_OUT );
750         if ( isIN && isOUT )
751           (*ip1)._transition = Trans_TANGENT;
752         else
753           (*ip1)._transition = isIN ? Trans_IN : Trans_OUT;
754       }
755     }
756   }
757   //================================================================================
758   /*
759    * Return "is OUT" state for nodes before the given intersection point
760    */
761   bool GridLine::GetIsOutBefore( multiset< F_IntersectPoint >::iterator ip, bool prevIsOut )
762   {
763     if ( ip->_transition == Trans_IN )
764       return true;
765     if ( ip->_transition == Trans_OUT )
766       return false;
767     if ( ip->_transition == Trans_APEX )
768     {
769       // singularity point (apex of a cone)
770       if ( _intPoints.size() == 1 || ip == _intPoints.begin() )
771         return true;
772       multiset< F_IntersectPoint >::iterator ipBef = ip, ipAft = ++ip;
773       if ( ipAft == _intPoints.end() )
774         return false;
775       --ipBef;
776       if ( ipBef->_transition != ipAft->_transition )
777         return ( ipBef->_transition == Trans_OUT );
778       return ( ipBef->_transition != Trans_OUT );
779     }
780     // _transition == Trans_TANGENT
781     return !prevIsOut;
782   }
783   //================================================================================
784   /*
785    * Adds face IDs
786    */
787   void B_IntersectPoint::Add( const vector< TGeomID >& fIDs,
788                               const SMDS_MeshNode*     n) const
789   {
790     if ( _faceIDs.empty() )
791       _faceIDs = fIDs;
792     else
793       for ( size_t i = 0; i < fIDs.size(); ++i )
794       {
795         vector< TGeomID >::iterator it =
796           std::find( _faceIDs.begin(), _faceIDs.end(), fIDs[i] );
797         if ( it == _faceIDs.end() )
798           _faceIDs.push_back( fIDs[i] );
799       }
800     if ( !_node )
801       _node = n;
802   }
803   //================================================================================
804   /*
805    * Returns index of a common face if any, else zero
806    */
807   int B_IntersectPoint::HasCommonFace( const B_IntersectPoint * other, int avoidFace ) const
808   {
809     if ( other )
810       for ( size_t i = 0; i < other->_faceIDs.size(); ++i )
811         if ( avoidFace != other->_faceIDs[i] &&
812              IsOnFace   ( other->_faceIDs[i] ))
813           return other->_faceIDs[i];
814     return 0;
815   }
816   //================================================================================
817   /*
818    * Returns \c true if \a faceID in in this->_faceIDs
819    */
820   bool B_IntersectPoint::IsOnFace( int faceID ) const // returns true if faceID is found
821   {
822     vector< TGeomID >::const_iterator it =
823       std::find( _faceIDs.begin(), _faceIDs.end(), faceID );
824     return ( it != _faceIDs.end() );
825   }
826   //================================================================================
827   /*
828    * Return an iterator on GridLine's in a given direction
829    */
830   LineIndexer Grid::GetLineIndexer(size_t iDir) const
831   {
832     const size_t indices[] = { 1,2,0, 0,2,1, 0,1,2 };
833     const string s      [] = { "X", "Y", "Z" };
834     LineIndexer li( _coords[0].size(),  _coords[1].size(),    _coords[2].size(),
835                     indices[iDir*3],    indices[iDir*3+1],    indices[iDir*3+2],
836                     s[indices[iDir*3]], s[indices[iDir*3+1]], s[indices[iDir*3+2]]);
837     return li;
838   }
839   //=============================================================================
840   /*
841    * Creates GridLine's of the grid
842    */
843   void Grid::SetCoordinates(const vector<double>& xCoords,
844                             const vector<double>& yCoords,
845                             const vector<double>& zCoords,
846                             const double*         axesDirs,
847                             const Bnd_Box&        shapeBox)
848   {
849     _coords[0] = xCoords;
850     _coords[1] = yCoords;
851     _coords[2] = zCoords;
852
853     _axes[0].SetCoord( axesDirs[0],
854                        axesDirs[1],
855                        axesDirs[2]);
856     _axes[1].SetCoord( axesDirs[3],
857                        axesDirs[4],
858                        axesDirs[5]);
859     _axes[2].SetCoord( axesDirs[6],
860                        axesDirs[7],
861                        axesDirs[8]);
862     _axes[0].Normalize();
863     _axes[1].Normalize();
864     _axes[2].Normalize();
865
866     _invB.SetCols( _axes[0], _axes[1], _axes[2] );
867     _invB.Invert();
868
869     // compute tolerance
870     _minCellSize = Precision::Infinite();
871     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
872     {
873       for ( size_t i = 1; i < _coords[ iDir ].size(); ++i )
874       {
875         double cellLen = _coords[ iDir ][ i ] - _coords[ iDir ][ i-1 ];
876         if ( cellLen < _minCellSize )
877           _minCellSize = cellLen;
878       }
879     }
880     if ( _minCellSize < Precision::Confusion() )
881       throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
882                                 SMESH_Comment("Too small cell size: ") << _minCellSize );
883     _tol = _minCellSize / 1000.;
884
885     // attune grid extremities to shape bounding box
886
887     double sP[6]; // aXmin, aYmin, aZmin, aXmax, aYmax, aZmax
888     shapeBox.Get(sP[0],sP[1],sP[2],sP[3],sP[4],sP[5]);
889     double* cP[6] = { &_coords[0].front(), &_coords[1].front(), &_coords[2].front(),
890                       &_coords[0].back(),  &_coords[1].back(),  &_coords[2].back() };
891     for ( int i = 0; i < 6; ++i )
892       if ( fabs( sP[i] - *cP[i] ) < _tol )
893         *cP[i] = sP[i];// + _tol/1000. * ( i < 3 ? +1 : -1 );
894
895     for ( int iDir = 0; iDir < 3; ++iDir )
896     {
897       if ( _coords[iDir][0] - sP[iDir] > _tol )
898       {
899         _minCellSize = Min( _minCellSize, _coords[iDir][0] - sP[iDir] );
900         _coords[iDir].insert( _coords[iDir].begin(), sP[iDir] + _tol/1000.);
901       }
902       if ( sP[iDir+3] - _coords[iDir].back() > _tol  )
903       {
904         _minCellSize = Min( _minCellSize, sP[iDir+3] - _coords[iDir].back() );
905         _coords[iDir].push_back( sP[iDir+3] - _tol/1000.);
906       }
907     }
908     _tol = _minCellSize / 1000.;
909
910     _origin = ( _coords[0][0] * _axes[0] +
911                 _coords[1][0] * _axes[1] +
912                 _coords[2][0] * _axes[2] );
913
914     // create lines
915     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
916     {
917       LineIndexer li = GetLineIndexer( iDir );
918       _lines[iDir].resize( li.NbLines() );
919       double len = _coords[ iDir ].back() - _coords[iDir].front();
920       for ( ; li.More(); ++li )
921       {
922         GridLine& gl = _lines[iDir][ li.LineIndex() ];
923         gl._line.SetLocation( _coords[0][li.I()] * _axes[0] +
924                               _coords[1][li.J()] * _axes[1] +
925                               _coords[2][li.K()] * _axes[2] );
926         gl._line.SetDirection( _axes[ iDir ]);
927         gl._length = len;
928       }
929     }
930   }
931   //================================================================================
932   /*
933    * Computes coordinates of a point in the grid CS
934    */
935   void Grid::ComputeUVW(const gp_XYZ& P, double UVW[3])
936   {
937     gp_XYZ p = P * _invB;
938     p.Coord( UVW[0], UVW[1], UVW[2] );
939   }
940   //================================================================================
941   /*
942    * Creates all nodes
943    */
944   void Grid::ComputeNodes(SMESH_MesherHelper& helper)
945   {
946     // state of each node of the grid relative to the geometry
947     const size_t nbGridNodes = _coords[0].size() * _coords[1].size() * _coords[2].size();
948     vector< bool > isNodeOut( nbGridNodes, false );
949     _nodes.resize( nbGridNodes, 0 );
950     _gridIntP.resize( nbGridNodes, NULL );
951
952     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
953     {
954       LineIndexer li = GetLineIndexer( iDir );
955
956       // find out a shift of node index while walking along a GridLine in this direction
957       li.SetIndexOnLine( 0 );
958       size_t nIndex0 = NodeIndex( li.I(), li.J(), li.K() );
959       li.SetIndexOnLine( 1 );
960       const size_t nShift = NodeIndex( li.I(), li.J(), li.K() ) - nIndex0;
961       
962       const vector<double> & coords = _coords[ iDir ];
963       for ( ; li.More(); ++li ) // loop on lines in iDir
964       {
965         li.SetIndexOnLine( 0 );
966         nIndex0 = NodeIndex( li.I(), li.J(), li.K() );
967
968         GridLine& line = _lines[ iDir ][ li.LineIndex() ];
969         const gp_XYZ lineLoc = line._line.Location().XYZ();
970         const gp_XYZ lineDir = line._line.Direction().XYZ();
971         line.RemoveExcessIntPoints( _tol );
972         multiset< F_IntersectPoint >& intPnts = line._intPoints;
973         multiset< F_IntersectPoint >::iterator ip = intPnts.begin();
974
975         bool isOut = true;
976         const double* nodeCoord = & coords[0];
977         const double* coord0    = nodeCoord;
978         const double* coordEnd  = coord0 + coords.size();
979         double nodeParam = 0;
980         for ( ; ip != intPnts.end(); ++ip )
981         {
982           // set OUT state or just skip IN nodes before ip
983           if ( nodeParam < ip->_paramOnLine - _tol )
984           {
985             isOut = line.GetIsOutBefore( ip, isOut );
986
987             while ( nodeParam < ip->_paramOnLine - _tol )
988             {
989               if ( isOut )
990                 isNodeOut[ nIndex0 + nShift * ( nodeCoord-coord0 ) ] = isOut;
991               if ( ++nodeCoord <  coordEnd )
992                 nodeParam = *nodeCoord - *coord0;
993               else
994                 break;
995             }
996             if ( nodeCoord == coordEnd ) break;
997           }
998           // create a mesh node on a GridLine at ip if it does not coincide with a grid node
999           if ( nodeParam > ip->_paramOnLine + _tol )
1000           {
1001             // li.SetIndexOnLine( 0 );
1002             // double xyz[3] = { _coords[0][ li.I() ], _coords[1][ li.J() ], _coords[2][ li.K() ]};
1003             // xyz[ li._iConst ] += ip->_paramOnLine;
1004             gp_XYZ xyz = lineLoc + ip->_paramOnLine * lineDir;
1005             ip->_node = helper.AddNode( xyz.X(), xyz.Y(), xyz.Z() );
1006             ip->_indexOnLine = nodeCoord-coord0-1;
1007           }
1008           // create a mesh node at ip concident with a grid node
1009           else
1010           {
1011             int nodeIndex = nIndex0 + nShift * ( nodeCoord-coord0 );
1012             if ( !_nodes[ nodeIndex ] )
1013             {
1014               //li.SetIndexOnLine( nodeCoord-coord0 );
1015               //double xyz[3] = { _coords[0][ li.I() ], _coords[1][ li.J() ], _coords[2][ li.K() ]};
1016               gp_XYZ xyz = lineLoc + nodeParam * lineDir;
1017               _nodes   [ nodeIndex ] = helper.AddNode( xyz.X(), xyz.Y(), xyz.Z() );
1018               _gridIntP[ nodeIndex ] = & * ip;
1019             }
1020             if ( _gridIntP[ nodeIndex ] )
1021               _gridIntP[ nodeIndex ]->Add( ip->_faceIDs );
1022             else
1023               _gridIntP[ nodeIndex ] = & * ip;
1024             // ip->_node        = _nodes[ nodeIndex ]; -- to differ from ip on links
1025             ip->_indexOnLine = nodeCoord-coord0;
1026             if ( ++nodeCoord < coordEnd )
1027               nodeParam = *nodeCoord - *coord0;
1028           }
1029         }
1030         // set OUT state to nodes after the last ip
1031         for ( ; nodeCoord < coordEnd; ++nodeCoord )
1032           isNodeOut[ nIndex0 + nShift * ( nodeCoord-coord0 ) ] = true;
1033       }
1034     }
1035
1036     // Create mesh nodes at !OUT nodes of the grid
1037
1038     for ( size_t z = 0; z < _coords[2].size(); ++z )
1039       for ( size_t y = 0; y < _coords[1].size(); ++y )
1040         for ( size_t x = 0; x < _coords[0].size(); ++x )
1041         {
1042           size_t nodeIndex = NodeIndex( x, y, z );
1043           if ( !isNodeOut[ nodeIndex ] && !_nodes[ nodeIndex] )
1044           {
1045             //_nodes[ nodeIndex ] = helper.AddNode( _coords[0][x], _coords[1][y], _coords[2][z] );
1046             gp_XYZ xyz = ( _coords[0][x] * _axes[0] +
1047                            _coords[1][y] * _axes[1] +
1048                            _coords[2][z] * _axes[2] );
1049             _nodes[ nodeIndex ] = helper.AddNode( xyz.X(), xyz.Y(), xyz.Z() );
1050           }
1051         }
1052
1053 #ifdef _MY_DEBUG_
1054     // check validity of transitions
1055     const char* trName[] = { "TANGENT", "IN", "OUT", "APEX" };
1056     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
1057     {
1058       LineIndexer li = GetLineIndexer( iDir );
1059       for ( ; li.More(); ++li )
1060       {
1061         multiset< F_IntersectPoint >& intPnts = _lines[ iDir ][ li.LineIndex() ]._intPoints;
1062         if ( intPnts.empty() ) continue;
1063         if ( intPnts.size() == 1 )
1064         {
1065           if ( intPnts.begin()->_transition != Trans_TANGENT &&
1066                intPnts.begin()->_transition != Trans_APEX )
1067           throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
1068                                     SMESH_Comment("Wrong SOLE transition of GridLine (")
1069                                     << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
1070                                     << ") along " << li._nameConst
1071                                     << ": " << trName[ intPnts.begin()->_transition] );
1072         }
1073         else
1074         {
1075           if ( intPnts.begin()->_transition == Trans_OUT )
1076             throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
1077                                       SMESH_Comment("Wrong START transition of GridLine (")
1078                                       << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
1079                                       << ") along " << li._nameConst
1080                                       << ": " << trName[ intPnts.begin()->_transition ]);
1081           if ( intPnts.rbegin()->_transition == Trans_IN )
1082             throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
1083                                       SMESH_Comment("Wrong END transition of GridLine (")
1084                                       << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
1085                                       << ") along " << li._nameConst
1086                                     << ": " << trName[ intPnts.rbegin()->_transition ]);
1087         }
1088       }
1089     }
1090 #endif
1091   }
1092
1093   //=============================================================================
1094   /*
1095    * Intersects TopoDS_Face with all GridLine's
1096    */
1097   void FaceGridIntersector::Intersect()
1098   {
1099     FaceLineIntersector intersector;
1100     intersector._surfaceInt = GetCurveFaceIntersector();
1101     intersector._tol        = _grid->_tol;
1102     intersector._transOut   = _face.Orientation() == TopAbs_REVERSED ? Trans_IN : Trans_OUT;
1103     intersector._transIn    = _face.Orientation() == TopAbs_REVERSED ? Trans_OUT : Trans_IN;
1104
1105     typedef void (FaceLineIntersector::* PIntFun )(const GridLine& gridLine);
1106     PIntFun interFunction;
1107
1108     bool isDirect = true;
1109     BRepAdaptor_Surface surf( _face );
1110     switch ( surf.GetType() ) {
1111     case GeomAbs_Plane:
1112       intersector._plane = surf.Plane();
1113       interFunction = &FaceLineIntersector::IntersectWithPlane;
1114       isDirect = intersector._plane.Direct();
1115       break;
1116     case GeomAbs_Cylinder:
1117       intersector._cylinder = surf.Cylinder();
1118       interFunction = &FaceLineIntersector::IntersectWithCylinder;
1119       isDirect = intersector._cylinder.Direct();
1120       break;
1121     case GeomAbs_Cone:
1122       intersector._cone = surf.Cone();
1123       interFunction = &FaceLineIntersector::IntersectWithCone;
1124       //isDirect = intersector._cone.Direct();
1125       break;
1126     case GeomAbs_Sphere:
1127       intersector._sphere = surf.Sphere();
1128       interFunction = &FaceLineIntersector::IntersectWithSphere;
1129       isDirect = intersector._sphere.Direct();
1130       break;
1131     case GeomAbs_Torus:
1132       intersector._torus = surf.Torus();
1133       interFunction = &FaceLineIntersector::IntersectWithTorus;
1134       //isDirect = intersector._torus.Direct();
1135       break;
1136     default:
1137       interFunction = &FaceLineIntersector::IntersectWithSurface;
1138     }
1139     if ( !isDirect )
1140       std::swap( intersector._transOut, intersector._transIn );
1141
1142     _intersections.clear();
1143     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
1144     {
1145       if ( surf.GetType() == GeomAbs_Plane )
1146       {
1147         // check if all lines in this direction are parallel to a plane
1148         if ( intersector._plane.Axis().IsNormal( _grid->_lines[iDir][0]._line.Position(),
1149                                                  Precision::Angular()))
1150           continue;
1151         // find out a transition, that is the same for all lines of a direction
1152         gp_Dir plnNorm = intersector._plane.Axis().Direction();
1153         gp_Dir lineDir = _grid->_lines[iDir][0]._line.Direction();
1154         intersector._transition =
1155           ( plnNorm * lineDir < 0 ) ? intersector._transIn : intersector._transOut;
1156       }
1157       if ( surf.GetType() == GeomAbs_Cylinder )
1158       {
1159         // check if all lines in this direction are parallel to a cylinder
1160         if ( intersector._cylinder.Axis().IsParallel( _grid->_lines[iDir][0]._line.Position(),
1161                                                       Precision::Angular()))
1162           continue;
1163       }
1164
1165       // intersect the grid lines with the face
1166       for ( size_t iL = 0; iL < _grid->_lines[iDir].size(); ++iL )
1167       {
1168         GridLine& gridLine = _grid->_lines[iDir][iL];
1169         if ( _bndBox.IsOut( gridLine._line )) continue;
1170
1171         intersector._intPoints.clear();
1172         (intersector.*interFunction)( gridLine ); // <- intersection with gridLine
1173         for ( size_t i = 0; i < intersector._intPoints.size(); ++i )
1174           _intersections.push_back( make_pair( &gridLine, intersector._intPoints[i] ));
1175       }
1176     }
1177   }
1178   //================================================================================
1179   /*
1180    * Return true if (_u,_v) is on the face
1181    */
1182   bool FaceLineIntersector::UVIsOnFace() const
1183   {
1184     TopAbs_State state = _surfaceInt->ClassifyUVPoint(gp_Pnt2d( _u,_v ));
1185     return ( state == TopAbs_IN || state == TopAbs_ON );
1186   }
1187   //================================================================================
1188   /*
1189    * Store an intersection if it is IN or ON the face
1190    */
1191   void FaceLineIntersector::addIntPoint(const bool toClassify)
1192   {
1193     if ( !toClassify || UVIsOnFace() )
1194     {
1195       F_IntersectPoint p;
1196       p._paramOnLine = _w;
1197       p._transition  = _transition;
1198       _intPoints.push_back( p );
1199     }
1200   }
1201   //================================================================================
1202   /*
1203    * Intersect a line with a plane
1204    */
1205   void FaceLineIntersector::IntersectWithPlane(const GridLine& gridLine)
1206   {
1207     IntAna_IntConicQuad linPlane( gridLine._line, _plane, Precision::Angular());
1208     _w = linPlane.ParamOnConic(1);
1209     if ( isParamOnLineOK( gridLine._length ))
1210     {
1211       ElSLib::Parameters(_plane, linPlane.Point(1) ,_u,_v);
1212       addIntPoint();
1213     }
1214   }
1215   //================================================================================
1216   /*
1217    * Intersect a line with a cylinder
1218    */
1219   void FaceLineIntersector::IntersectWithCylinder(const GridLine& gridLine)
1220   {
1221     IntAna_IntConicQuad linCylinder( gridLine._line, _cylinder );
1222     if ( linCylinder.IsDone() && linCylinder.NbPoints() > 0 )
1223     {
1224       _w = linCylinder.ParamOnConic(1);
1225       if ( linCylinder.NbPoints() == 1 )
1226         _transition = Trans_TANGENT;
1227       else
1228         _transition = _w < linCylinder.ParamOnConic(2) ? _transIn : _transOut;
1229       if ( isParamOnLineOK( gridLine._length ))
1230       {
1231         ElSLib::Parameters(_cylinder, linCylinder.Point(1) ,_u,_v);
1232         addIntPoint();
1233       }
1234       if ( linCylinder.NbPoints() > 1 )
1235       {
1236         _w = linCylinder.ParamOnConic(2);
1237         if ( isParamOnLineOK( gridLine._length ))
1238         {
1239           ElSLib::Parameters(_cylinder, linCylinder.Point(2) ,_u,_v);
1240           _transition = ( _transition == Trans_OUT ) ? Trans_IN : Trans_OUT;
1241           addIntPoint();
1242         }
1243       }
1244     }
1245   }
1246   //================================================================================
1247   /*
1248    * Intersect a line with a cone
1249    */
1250   void FaceLineIntersector::IntersectWithCone (const GridLine& gridLine)
1251   {
1252     IntAna_IntConicQuad linCone(gridLine._line,_cone);
1253     if ( !linCone.IsDone() ) return;
1254     gp_Pnt P;
1255     gp_Vec du, dv, norm;
1256     for ( int i = 1; i <= linCone.NbPoints(); ++i )
1257     {
1258       _w = linCone.ParamOnConic( i );
1259       if ( !isParamOnLineOK( gridLine._length )) continue;
1260       ElSLib::Parameters(_cone, linCone.Point(i) ,_u,_v);
1261       if ( UVIsOnFace() )
1262       {
1263         ElSLib::D1( _u, _v, _cone, P, du, dv );
1264         norm = du ^ dv;
1265         double normSize2 = norm.SquareMagnitude();
1266         if ( normSize2 > Precision::Angular() * Precision::Angular() )
1267         {
1268           double cos = norm.XYZ() * gridLine._line.Direction().XYZ();
1269           cos /= sqrt( normSize2 );
1270           if ( cos < -Precision::Angular() )
1271             _transition = _transIn;
1272           else if ( cos > Precision::Angular() )
1273             _transition = _transOut;
1274           else
1275             _transition = Trans_TANGENT;
1276         }
1277         else
1278         {
1279           _transition = Trans_APEX;
1280         }
1281         addIntPoint( /*toClassify=*/false);
1282       }
1283     }
1284   }
1285   //================================================================================
1286   /*
1287    * Intersect a line with a sphere
1288    */
1289   void FaceLineIntersector::IntersectWithSphere  (const GridLine& gridLine)
1290   {
1291     IntAna_IntConicQuad linSphere(gridLine._line,_sphere);
1292     if ( linSphere.IsDone() && linSphere.NbPoints() > 0 )
1293     {
1294       _w = linSphere.ParamOnConic(1);
1295       if ( linSphere.NbPoints() == 1 )
1296         _transition = Trans_TANGENT;
1297       else
1298         _transition = _w < linSphere.ParamOnConic(2) ? _transIn : _transOut;
1299       if ( isParamOnLineOK( gridLine._length ))
1300       {
1301         ElSLib::Parameters(_sphere, linSphere.Point(1) ,_u,_v);
1302         addIntPoint();
1303       }
1304       if ( linSphere.NbPoints() > 1 )
1305       {
1306         _w = linSphere.ParamOnConic(2);
1307         if ( isParamOnLineOK( gridLine._length ))
1308         {
1309           ElSLib::Parameters(_sphere, linSphere.Point(2) ,_u,_v);
1310           _transition = ( _transition == Trans_OUT ) ? Trans_IN : Trans_OUT;
1311           addIntPoint();
1312         }
1313       }
1314     }
1315   }
1316   //================================================================================
1317   /*
1318    * Intersect a line with a torus
1319    */
1320   void FaceLineIntersector::IntersectWithTorus   (const GridLine& gridLine)
1321   {
1322     IntAna_IntLinTorus linTorus(gridLine._line,_torus);
1323     if ( !linTorus.IsDone()) return;
1324     gp_Pnt P;
1325     gp_Vec du, dv, norm;
1326     for ( int i = 1; i <= linTorus.NbPoints(); ++i )
1327     {
1328       _w = linTorus.ParamOnLine( i );
1329       if ( !isParamOnLineOK( gridLine._length )) continue;
1330       linTorus.ParamOnTorus( i, _u,_v );
1331       if ( UVIsOnFace() )
1332       {
1333         ElSLib::D1( _u, _v, _torus, P, du, dv );
1334         norm = du ^ dv;
1335         double normSize = norm.Magnitude();
1336         double cos = norm.XYZ() * gridLine._line.Direction().XYZ();
1337         cos /= normSize;
1338         if ( cos < -Precision::Angular() )
1339           _transition = _transIn;
1340         else if ( cos > Precision::Angular() )
1341           _transition = _transOut;
1342         else
1343           _transition = Trans_TANGENT;
1344         addIntPoint( /*toClassify=*/false);
1345       }
1346     }
1347   }
1348   //================================================================================
1349   /*
1350    * Intersect a line with a non-analytical surface
1351    */
1352   void FaceLineIntersector::IntersectWithSurface (const GridLine& gridLine)
1353   {
1354     _surfaceInt->Perform( gridLine._line, 0.0, gridLine._length );
1355     if ( !_surfaceInt->IsDone() ) return;
1356     for ( int i = 1; i <= _surfaceInt->NbPnt(); ++i )
1357     {
1358       _transition = Transition( _surfaceInt->Transition( i ) );
1359       _w = _surfaceInt->WParameter( i );
1360       addIntPoint(/*toClassify=*/false);
1361     }
1362   }
1363   //================================================================================
1364   /*
1365    * check if its face can be safely intersected in a thread
1366    */
1367   bool FaceGridIntersector::IsThreadSafe(set< const Standard_Transient* >& noSafeTShapes) const
1368   {
1369     bool isSafe = true;
1370
1371     // check surface
1372     TopLoc_Location loc;
1373     Handle(Geom_Surface) surf = BRep_Tool::Surface( _face, loc );
1374     Handle(Geom_RectangularTrimmedSurface) ts =
1375       Handle(Geom_RectangularTrimmedSurface)::DownCast( surf );
1376     while( !ts.IsNull() ) {
1377       surf = ts->BasisSurface();
1378       ts = Handle(Geom_RectangularTrimmedSurface)::DownCast(surf);
1379     }
1380     if ( surf->IsKind( STANDARD_TYPE(Geom_BSplineSurface )) ||
1381          surf->IsKind( STANDARD_TYPE(Geom_BezierSurface )))
1382       if ( !noSafeTShapes.insert( _face.TShape().get() ).second )
1383         isSafe = false;
1384
1385     double f, l;
1386     TopExp_Explorer exp( _face, TopAbs_EDGE );
1387     for ( ; exp.More(); exp.Next() )
1388     {
1389       bool edgeIsSafe = true;
1390       const TopoDS_Edge& e = TopoDS::Edge( exp.Current() );
1391       // check 3d curve
1392       {
1393         Handle(Geom_Curve) c = BRep_Tool::Curve( e, loc, f, l);
1394         if ( !c.IsNull() )
1395         {
1396           Handle(Geom_TrimmedCurve) tc = Handle(Geom_TrimmedCurve)::DownCast(c);
1397           while( !tc.IsNull() ) {
1398             c = tc->BasisCurve();
1399             tc = Handle(Geom_TrimmedCurve)::DownCast(c);
1400           }
1401           if ( c->IsKind( STANDARD_TYPE(Geom_BSplineCurve )) ||
1402                c->IsKind( STANDARD_TYPE(Geom_BezierCurve )))
1403             edgeIsSafe = false;
1404         }
1405       }
1406       // check 2d curve
1407       if ( edgeIsSafe )
1408       {
1409         Handle(Geom2d_Curve) c2 = BRep_Tool::CurveOnSurface( e, surf, loc, f, l);
1410         if ( !c2.IsNull() )
1411         {
1412           Handle(Geom2d_TrimmedCurve) tc = Handle(Geom2d_TrimmedCurve)::DownCast(c2);
1413           while( !tc.IsNull() ) {
1414             c2 = tc->BasisCurve();
1415             tc = Handle(Geom2d_TrimmedCurve)::DownCast(c2);
1416           }
1417           if ( c2->IsKind( STANDARD_TYPE(Geom2d_BSplineCurve )) ||
1418                c2->IsKind( STANDARD_TYPE(Geom2d_BezierCurve )))
1419             edgeIsSafe = false;
1420         }
1421       }
1422       if ( !edgeIsSafe && !noSafeTShapes.insert( e.TShape().get() ).second )
1423         isSafe = false;
1424     }
1425     return isSafe;
1426   }
1427   //================================================================================
1428   /*!
1429    * \brief Creates topology of the hexahedron
1430    */
1431   Hexahedron::Hexahedron(const double sizeThreshold, Grid* grid)
1432     : _grid( grid ), _sizeThreshold( sizeThreshold ), _nbFaceIntNodes(0)
1433   {
1434     _polygons.reserve(100); // to avoid reallocation;
1435
1436     //set nodes shift within grid->_nodes from the node 000 
1437     size_t dx = _grid->NodeIndexDX();
1438     size_t dy = _grid->NodeIndexDY();
1439     size_t dz = _grid->NodeIndexDZ();
1440     size_t i000 = 0;
1441     size_t i100 = i000 + dx;
1442     size_t i010 = i000 + dy;
1443     size_t i110 = i010 + dx;
1444     size_t i001 = i000 + dz;
1445     size_t i101 = i100 + dz;
1446     size_t i011 = i010 + dz;
1447     size_t i111 = i110 + dz;
1448     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V000 )] = i000;
1449     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V100 )] = i100;
1450     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V010 )] = i010;
1451     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V110 )] = i110;
1452     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V001 )] = i001;
1453     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V101 )] = i101;
1454     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V011 )] = i011;
1455     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V111 )] = i111;
1456
1457     vector< int > idVec;
1458     // set nodes to links
1459     for ( int linkID = SMESH_Block::ID_Ex00; linkID <= SMESH_Block::ID_E11z; ++linkID )
1460     {
1461       SMESH_Block::GetEdgeVertexIDs( linkID, idVec );
1462       _Link& link = _hexLinks[ SMESH_Block::ShapeIndex( linkID )];
1463       link._nodes[0] = &_hexNodes[ SMESH_Block::ShapeIndex( idVec[0] )];
1464       link._nodes[1] = &_hexNodes[ SMESH_Block::ShapeIndex( idVec[1] )];
1465     }
1466
1467     // set links to faces
1468     int interlace[4] = { 0, 3, 1, 2 }; // to walk by links around a face: { u0, 1v, u1, 0v }
1469     for ( int faceID = SMESH_Block::ID_Fxy0; faceID <= SMESH_Block::ID_F1yz; ++faceID )
1470     {
1471       SMESH_Block::GetFaceEdgesIDs( faceID, idVec );
1472       _Face& quad = _hexQuads[ SMESH_Block::ShapeIndex( faceID )];
1473       bool revFace = ( faceID == SMESH_Block::ID_Fxy0 ||
1474                        faceID == SMESH_Block::ID_Fx1z ||
1475                        faceID == SMESH_Block::ID_F0yz );
1476       quad._links.resize(4);
1477       vector<_OrientedLink>::iterator         frwLinkIt = quad._links.begin();
1478       vector<_OrientedLink>::reverse_iterator revLinkIt = quad._links.rbegin();
1479       for ( int i = 0; i < 4; ++i )
1480       {
1481         bool revLink = revFace;
1482         if ( i > 1 ) // reverse links u1 and v0
1483           revLink = !revLink;
1484         _OrientedLink& link = revFace ? *revLinkIt++ : *frwLinkIt++;
1485         link = _OrientedLink( & _hexLinks[ SMESH_Block::ShapeIndex( idVec[interlace[i]] )],
1486                               revLink );
1487       }
1488     }
1489   }
1490   //================================================================================
1491   /*!
1492    * \brief Copy constructor
1493    */
1494   Hexahedron::Hexahedron( const Hexahedron& other )
1495     :_grid( other._grid ), _sizeThreshold( other._sizeThreshold ), _nbFaceIntNodes(0)
1496   {
1497     _polygons.reserve(100); // to avoid reallocation;
1498
1499     for ( int i = 0; i < 8; ++i )
1500       _nodeShift[i] = other._nodeShift[i];
1501
1502     for ( int i = 0; i < 12; ++i )
1503     {
1504       const _Link& srcLink = other._hexLinks[ i ];
1505       _Link&       tgtLink = this->_hexLinks[ i ];
1506       tgtLink._nodes[0] = _hexNodes + ( srcLink._nodes[0] - other._hexNodes );
1507       tgtLink._nodes[1] = _hexNodes + ( srcLink._nodes[1] - other._hexNodes );
1508     }
1509
1510     for ( int i = 0; i < 6; ++i )
1511     {
1512       const _Face& srcQuad = other._hexQuads[ i ];
1513       _Face&       tgtQuad = this->_hexQuads[ i ];
1514       tgtQuad._links.resize(4);
1515       for ( int j = 0; j < 4; ++j )
1516       {
1517         const _OrientedLink& srcLink = srcQuad._links[ j ];
1518         _OrientedLink&       tgtLink = tgtQuad._links[ j ];
1519         tgtLink._reverse = srcLink._reverse;
1520         tgtLink._link    = _hexLinks + ( srcLink._link - other._hexLinks );
1521       }
1522     }
1523   }
1524
1525   //================================================================================
1526   /*!
1527    * \brief Initializes its data by given grid cell
1528    */
1529   void Hexahedron::init( size_t i, size_t j, size_t k )
1530   {
1531     _i = i; _j = j; _k = k;
1532     // set nodes of grid to nodes of the hexahedron and
1533     // count nodes at hexahedron corners located IN and ON geometry
1534     _nbCornerNodes = _nbBndNodes = 0;
1535     _origNodeInd   = _grid->NodeIndex( i,j,k );
1536     for ( int iN = 0; iN < 8; ++iN )
1537     {
1538       _hexNodes[iN]._node     = _grid->_nodes   [ _origNodeInd + _nodeShift[iN] ];
1539       _hexNodes[iN]._intPoint = _grid->_gridIntP[ _origNodeInd + _nodeShift[iN] ];
1540       _nbCornerNodes += bool( _hexNodes[iN]._node );
1541       _nbBndNodes    += bool( _hexNodes[iN]._intPoint );
1542     }
1543     _sideLength[0] = _grid->_coords[0][i+1] - _grid->_coords[0][i];
1544     _sideLength[1] = _grid->_coords[1][j+1] - _grid->_coords[1][j];
1545     _sideLength[2] = _grid->_coords[2][k+1] - _grid->_coords[2][k];
1546
1547     _intNodes.clear();
1548     _vIntNodes.clear();
1549
1550     if ( _nbFaceIntNodes + _eIntPoints.size() > 0 &&
1551          _nbFaceIntNodes + _nbCornerNodes + _eIntPoints.size() > 3)
1552     {
1553       _intNodes.reserve( 3 * _nbBndNodes + _nbFaceIntNodes + _eIntPoints.size() );
1554
1555       // this method can be called in parallel, so use own helper
1556       SMESH_MesherHelper helper( *_grid->_helper->GetMesh() );
1557
1558       // create sub-links (_splits) by splitting links with _fIntPoints
1559       _Link split;
1560       for ( int iLink = 0; iLink < 12; ++iLink )
1561       {
1562         _Link& link = _hexLinks[ iLink ];
1563         link._fIntNodes.resize( link._fIntPoints.size() );
1564         for ( size_t i = 0; i < link._fIntPoints.size(); ++i )
1565         {
1566           _intNodes.push_back( _Node( 0, link._fIntPoints[i] ));
1567           link._fIntNodes[ i ] = & _intNodes.back();
1568         }
1569
1570         link._splits.clear();
1571         split._nodes[ 0 ] = link._nodes[0];
1572         bool isOut = ( ! link._nodes[0]->Node() );
1573         bool checkTransition;
1574         for ( size_t i = 0; i < link._fIntNodes.size(); ++i )
1575         {
1576           const bool isGridNode = ( ! link._fIntNodes[i]->Node() );
1577           if ( !isGridNode ) // intersection non-coincident with a grid node
1578           {
1579             if ( split._nodes[ 0 ]->Node() && !isOut )
1580             {
1581               split._nodes[ 1 ] = link._fIntNodes[i];
1582               link._splits.push_back( split );
1583             }
1584             split._nodes[ 0 ] = link._fIntNodes[i];
1585             checkTransition = true;
1586           }
1587           else // FACE intersection coincident with a grid node (at link ends)
1588           {
1589             checkTransition = ( i == 0 && link._nodes[0]->Node() );
1590           }
1591           if ( checkTransition )
1592           {
1593             if ( link._fIntPoints[i]->_faceIDs.size() > 1 || _eIntPoints.size() > 0 )
1594               isOut = isOutPoint( link, i, helper );
1595             else
1596               switch ( link._fIntPoints[i]->_transition ) {
1597               case Trans_OUT: isOut = true;  break;
1598               case Trans_IN : isOut = false; break;
1599               default:
1600                 isOut = isOutPoint( link, i, helper );
1601               }
1602           }
1603         }
1604         if ( link._nodes[ 1 ]->Node() && split._nodes[ 0 ]->Node() && !isOut )
1605         {
1606           split._nodes[ 1 ] = link._nodes[1];
1607           link._splits.push_back( split );
1608         }
1609       }
1610
1611       // Create _Node's at intersections with EDGEs.
1612
1613       const double tol2 = _grid->_tol * _grid->_tol;
1614       int facets[3], nbFacets, subEntity;
1615
1616       for ( size_t iP = 0; iP < _eIntPoints.size(); ++iP )
1617       {
1618         nbFacets = getEntity( _eIntPoints[iP], facets, subEntity );
1619         _Node* equalNode = 0;
1620         switch( nbFacets ) {
1621         case 1: // in a _Face
1622         {
1623           _Face& quad = _hexQuads[ facets[0] - SMESH_Block::ID_FirstF ];
1624           equalNode = findEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
1625           if ( equalNode ) {
1626             equalNode->Add( _eIntPoints[ iP ] );
1627           }
1628           else {
1629             _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
1630             quad._eIntNodes.push_back( & _intNodes.back() );
1631           }
1632           break;
1633         }
1634         case 2: // on a _Link
1635         {
1636           _Link& link = _hexLinks[ subEntity - SMESH_Block::ID_FirstE ];
1637           if ( link._splits.size() > 0 )
1638           {
1639             equalNode = findEqualNode( link._fIntNodes, _eIntPoints[ iP ], tol2 );
1640             if ( equalNode )
1641               equalNode->Add( _eIntPoints[ iP ] );
1642           }
1643           else
1644           {
1645             _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
1646             for ( int iF = 0; iF < 2; ++iF )
1647             {
1648               _Face& quad = _hexQuads[ facets[iF] - SMESH_Block::ID_FirstF ];
1649               equalNode = findEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
1650               if ( equalNode ) {
1651                 equalNode->Add( _eIntPoints[ iP ] );
1652               }
1653               else {
1654                 quad._eIntNodes.push_back( & _intNodes.back() );
1655               }
1656             }
1657           }
1658           break;
1659         }
1660         case 3: // at a corner
1661         {
1662           _Node& node = _hexNodes[ subEntity - SMESH_Block::ID_FirstV ];
1663           if ( node.Node() > 0 )
1664           {
1665             if ( node._intPoint )
1666               node._intPoint->Add( _eIntPoints[ iP ]->_faceIDs, _eIntPoints[ iP ]->_node );
1667           }
1668           else
1669           {
1670             _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
1671             for ( int iF = 0; iF < 3; ++iF )
1672             {
1673               _Face& quad = _hexQuads[ facets[iF] - SMESH_Block::ID_FirstF ];
1674               equalNode = findEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
1675               if ( equalNode ) {
1676                 equalNode->Add( _eIntPoints[ iP ] );
1677               }
1678               else {
1679                 quad._eIntNodes.push_back( & _intNodes.back() );
1680               }
1681             }
1682           }
1683           break;
1684         }
1685         } // switch( nbFacets )
1686
1687         if ( nbFacets == 0 ||
1688              _grid->_shapes( _eIntPoints[ iP ]->_shapeID ).ShapeType() == TopAbs_VERTEX )
1689         {
1690           equalNode = findEqualNode( _vIntNodes, _eIntPoints[ iP ], tol2 );
1691           if ( equalNode ) {
1692             equalNode->Add( _eIntPoints[ iP ] );
1693           }
1694           else if ( nbFacets == 0 ) {
1695             if ( _intNodes.empty() || _intNodes.back().EdgeIntPnt() != _eIntPoints[ iP ])
1696               _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
1697             _vIntNodes.push_back( & _intNodes.back() );
1698           }
1699         }
1700       } // loop on _eIntPoints
1701     }
1702     else if ( 3 < _nbCornerNodes && _nbCornerNodes < 8 ) // _nbFaceIntNodes == 0
1703     {
1704       _Link split;
1705       // create sub-links (_splits) of whole links
1706       for ( int iLink = 0; iLink < 12; ++iLink )
1707       {
1708         _Link& link = _hexLinks[ iLink ];
1709         link._splits.clear();
1710         if ( link._nodes[ 0 ]->Node() && link._nodes[ 1 ]->Node() )
1711         {
1712           split._nodes[ 0 ] = link._nodes[0];
1713           split._nodes[ 1 ] = link._nodes[1];
1714           link._splits.push_back( split );
1715         }
1716       }
1717     }
1718
1719   }
1720   //================================================================================
1721   /*!
1722    * \brief Initializes its data by given grid cell (countered from zero)
1723    */
1724   void Hexahedron::init( size_t iCell )
1725   {
1726     size_t iNbCell = _grid->_coords[0].size() - 1;
1727     size_t jNbCell = _grid->_coords[1].size() - 1;
1728     _i = iCell % iNbCell;
1729     _j = ( iCell % ( iNbCell * jNbCell )) / iNbCell;
1730     _k = iCell / iNbCell / jNbCell;
1731     init( _i, _j, _k );
1732   }
1733
1734   //================================================================================
1735   /*!
1736    * \brief Compute mesh volumes resulted from intersection of the Hexahedron
1737    */
1738   void Hexahedron::ComputeElements()
1739   {
1740     Init();
1741
1742     int nbIntersections = _nbFaceIntNodes + _eIntPoints.size();
1743     if ( _nbCornerNodes + nbIntersections < 4 )
1744       return;
1745
1746     if ( _nbBndNodes == _nbCornerNodes && nbIntersections == 0 && isInHole() )
1747       return;
1748
1749     _polygons.clear();
1750     _polygons.reserve( 20 );
1751
1752     // Create polygons from quadrangles
1753     // --------------------------------
1754
1755     vector< _OrientedLink > splits;
1756     vector<_Node*>          chainNodes;
1757     _Face*                  coplanarPolyg;
1758
1759     bool hasEdgeIntersections = !_eIntPoints.empty();
1760
1761     for ( int iF = 0; iF < 6; ++iF ) // loop on 6 sides of a hexahedron
1762     {
1763       _Face& quad = _hexQuads[ iF ] ;
1764
1765       _polygons.resize( _polygons.size() + 1 );
1766       _Face* polygon = &_polygons.back();
1767       polygon->_polyLinks.reserve( 20 );
1768
1769       splits.clear();
1770       for ( int iE = 0; iE < 4; ++iE ) // loop on 4 sides of a quadrangle
1771         for ( int iS = 0; iS < quad._links[ iE ].NbResultLinks(); ++iS )
1772           splits.push_back( quad._links[ iE ].ResultLink( iS ));
1773
1774       // add splits of links to a polygon and add _polyLinks to make
1775       // polygon's boundary closed
1776
1777       int nbSplits = splits.size();
1778       if (( nbSplits == 1 ) &&
1779           ( quad._eIntNodes.empty() ||
1780             splits[0].FirstNode()->IsLinked( splits[0].LastNode()->_intPoint )))
1781           //( quad._eIntNodes.empty() || _nbCornerNodes + nbIntersections > 6 ))
1782         nbSplits = 0;
1783
1784 #ifdef _DEBUG_
1785       for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
1786         if ( quad._eIntNodes[ iP ]->IsUsedInFace( polygon ))
1787           quad._eIntNodes[ iP ]->_usedInFace = 0;
1788 #endif
1789       size_t nbUsedEdgeNodes = 0;
1790       _Face* prevPolyg = 0; // polygon previously created from this quad
1791
1792       while ( nbSplits > 0 )
1793       {
1794         size_t iS = 0;
1795         while ( !splits[ iS ] )
1796           ++iS;
1797
1798         if ( !polygon->_links.empty() )
1799         {
1800           _polygons.resize( _polygons.size() + 1 );
1801           polygon = &_polygons.back();
1802           polygon->_polyLinks.reserve( 20 );
1803         }
1804         polygon->_links.push_back( splits[ iS ] );
1805         splits[ iS++ ]._link = 0;
1806         --nbSplits;
1807
1808         _Node* nFirst = polygon->_links.back().FirstNode();
1809         _Node *n1,*n2 = polygon->_links.back().LastNode();
1810         for ( ; nFirst != n2 && iS < splits.size(); ++iS )
1811         {
1812           _OrientedLink& split = splits[ iS ];
1813           if ( !split ) continue;
1814
1815           n1 = split.FirstNode();
1816           if ( n1 == n2 &&
1817                n1->_intPoint &&
1818                n1->_intPoint->_faceIDs.size() > 1 )
1819           {
1820             // n1 is at intersection with EDGE
1821             if ( findChainOnEdge( splits, polygon->_links.back(), split, iS, quad, chainNodes ))
1822             {
1823               for ( size_t i = 1; i < chainNodes.size(); ++i )
1824                 polygon->AddPolyLink( chainNodes[i-1], chainNodes[i], prevPolyg );
1825               prevPolyg = polygon;
1826               n2 = chainNodes.back();
1827               continue;
1828             }
1829           }
1830           else if ( n1 != n2 )
1831           {
1832             // try to connect to intersections with EDGEs
1833             if ( quad._eIntNodes.size() > nbUsedEdgeNodes  &&
1834                  findChain( n2, n1, quad, chainNodes ))
1835             {
1836               for ( size_t i = 1; i < chainNodes.size(); ++i )
1837               {
1838                 polygon->AddPolyLink( chainNodes[i-1], chainNodes[i] );
1839                 nbUsedEdgeNodes += ( chainNodes[i]->IsUsedInFace( polygon ));
1840               }
1841               if ( chainNodes.back() != n1 )
1842               {
1843                 n2 = chainNodes.back();
1844                 --iS;
1845                 continue;
1846               }
1847             }
1848             // try to connect to a split ending on the same FACE
1849             else
1850             {
1851               _OrientedLink foundSplit;
1852               for ( size_t i = iS; i < splits.size() && !foundSplit; ++i )
1853                 if (( foundSplit = splits[ i ]) &&
1854                     ( n2->IsLinked( foundSplit.FirstNode()->_intPoint )))
1855                 {
1856                   iS = i - 1;
1857                 }
1858                 else
1859                 {
1860                   foundSplit._link = 0;
1861                 }
1862               if ( foundSplit )
1863               {
1864                 if ( n2 != foundSplit.FirstNode() )
1865                 {
1866                   polygon->AddPolyLink( n2, foundSplit.FirstNode() );
1867                   n2 = foundSplit.FirstNode();
1868                 }
1869                 continue;
1870               }
1871               else
1872               {
1873                 if ( n2->IsLinked( nFirst->_intPoint ))
1874                   break;
1875                 polygon->AddPolyLink( n2, n1, prevPolyg );
1876               }
1877             }
1878           } // if ( n1 != n2 )
1879
1880           polygon->_links.push_back( split );
1881           split._link = 0;
1882           --nbSplits;
1883           n2 = polygon->_links.back().LastNode();
1884
1885         } // loop on splits
1886
1887         if ( nFirst != n2 ) // close a polygon
1888         {
1889           if ( !findChain( n2, nFirst, quad, chainNodes ))
1890           {
1891             if ( !closePolygon( polygon, chainNodes ))
1892               if ( !isImplementEdges() )
1893                 chainNodes.push_back( nFirst );
1894           }
1895           for ( size_t i = 1; i < chainNodes.size(); ++i )
1896           {
1897             polygon->AddPolyLink( chainNodes[i-1], chainNodes[i], prevPolyg );
1898             nbUsedEdgeNodes += bool( chainNodes[i]->IsUsedInFace( polygon ));
1899           }
1900         }
1901
1902         if ( polygon->_links.size() < 3 && nbSplits > 0 )
1903         {
1904           polygon->_polyLinks.clear();
1905           polygon->_links.clear();
1906         }
1907       } // while ( nbSplits > 0 )
1908
1909       if ( polygon->_links.size() < 3 )
1910       {
1911         _polygons.pop_back();
1912       }
1913     }  // loop on 6 hexahedron sides
1914
1915     // Create polygons closing holes in a polyhedron
1916     // ----------------------------------------------
1917
1918     // clear _usedInFace
1919     for ( size_t iN = 0; iN < _intNodes.size(); ++iN )
1920       _intNodes[ iN ]._usedInFace = 0;
1921
1922     // add polygons to their links and mark used nodes
1923     for ( size_t iP = 0; iP < _polygons.size(); ++iP )
1924     {
1925       _Face& polygon = _polygons[ iP ];
1926       for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
1927       {
1928         polygon._links[ iL ].AddFace( &polygon );
1929         polygon._links[ iL ].FirstNode()->_usedInFace = &polygon;
1930       }
1931     }
1932     // find free links
1933     vector< _OrientedLink* > freeLinks;
1934     freeLinks.reserve(20);
1935     for ( size_t iP = 0; iP < _polygons.size(); ++iP )
1936     {
1937       _Face& polygon = _polygons[ iP ];
1938       for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
1939         if ( polygon._links[ iL ].NbFaces() < 2 )
1940           freeLinks.push_back( & polygon._links[ iL ]);
1941     }
1942     int nbFreeLinks = freeLinks.size();
1943     if ( nbFreeLinks == 1 ) return;
1944
1945     // put not used intersection nodes to _vIntNodes
1946     int nbVertexNodes = 0; // nb not used vertex nodes
1947     {
1948       for ( size_t iN = 0; iN < _vIntNodes.size(); ++iN )
1949         nbVertexNodes += ( !_vIntNodes[ iN ]->IsUsedInFace() );
1950
1951       const double tol = 1e-3 * Min( Min( _sideLength[0], _sideLength[1] ), _sideLength[0] );
1952       for ( size_t iN = _nbFaceIntNodes; iN < _intNodes.size(); ++iN )
1953       {
1954         if ( _intNodes[ iN ].IsUsedInFace() ) continue;
1955         if ( dynamic_cast< const F_IntersectPoint* >( _intNodes[ iN ]._intPoint )) continue;
1956         _Node* equalNode =
1957           findEqualNode( _vIntNodes, _intNodes[ iN ].EdgeIntPnt(), tol*tol );
1958         if ( !equalNode )
1959         {
1960           _vIntNodes.push_back( &_intNodes[ iN ]);
1961           ++nbVertexNodes;
1962         }
1963       }
1964     }
1965
1966     set<TGeomID> usedFaceIDs;
1967     vector< TGeomID > faces;
1968     TGeomID curFace = 0;
1969     const size_t nbQuadPolygons = _polygons.size();
1970     E_IntersectPoint ipTmp;
1971
1972     // create polygons by making closed chains of free links
1973     size_t iPolygon = _polygons.size();
1974     while ( nbFreeLinks > 0 )
1975     {
1976       if ( iPolygon == _polygons.size() )
1977       {
1978         _polygons.resize( _polygons.size() + 1 );
1979         _polygons[ iPolygon ]._polyLinks.reserve( 20 );
1980         _polygons[ iPolygon ]._links.reserve( 20 );
1981       }
1982       _Face& polygon = _polygons[ iPolygon ];
1983
1984       _OrientedLink* curLink = 0;
1985       _Node*         curNode;
1986       if (( !hasEdgeIntersections ) ||
1987           ( nbFreeLinks < 4 && nbVertexNodes == 0 ))
1988       {
1989         // get a remaining link to start from
1990         for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
1991           if (( curLink = freeLinks[ iL ] ))
1992             freeLinks[ iL ] = 0;
1993         polygon._links.push_back( *curLink );
1994         --nbFreeLinks;
1995         do
1996         {
1997           // find all links connected to curLink
1998           curNode = curLink->FirstNode();
1999           curLink = 0;
2000           for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2001             if ( freeLinks[ iL ] && freeLinks[ iL ]->LastNode() == curNode )
2002             {
2003               curLink = freeLinks[ iL ];
2004               freeLinks[ iL ] = 0;
2005               --nbFreeLinks;
2006               polygon._links.push_back( *curLink );
2007             }
2008         } while ( curLink );
2009       }
2010       else // there are intersections with EDGEs
2011       {
2012         // get a remaining link to start from, one lying on minimal nb of FACEs
2013         {
2014           typedef pair< TGeomID, int > TFaceOfLink;
2015           TFaceOfLink faceOfLink( -1, -1 );
2016           TFaceOfLink facesOfLink[3] = { faceOfLink, faceOfLink, faceOfLink };
2017           for ( size_t iL = 0; iL < freeLinks.size(); ++iL )
2018             if ( freeLinks[ iL ] )
2019             {
2020               faces = freeLinks[ iL ]->GetNotUsedFace( usedFaceIDs );
2021               if ( faces.size() == 1 )
2022               {
2023                 faceOfLink = TFaceOfLink( faces[0], iL );
2024                 if ( !freeLinks[ iL ]->HasEdgeNodes() )
2025                   break;
2026                 facesOfLink[0] = faceOfLink;
2027               }
2028               else if ( facesOfLink[0].first < 0 )
2029               {
2030                 faceOfLink = TFaceOfLink(( faces.empty() ? -1 : faces[0]), iL );
2031                 facesOfLink[ 1 + faces.empty() ] = faceOfLink;
2032               }
2033             }
2034           for ( int i = 0; faceOfLink.first < 0 && i < 3; ++i )
2035             faceOfLink = facesOfLink[i];
2036
2037           if ( faceOfLink.first < 0 ) // all faces used
2038           {
2039             for ( size_t iL = 0; iL < freeLinks.size() && faceOfLink.first < 1; ++iL )
2040               if (( curLink = freeLinks[ iL ]))
2041               {
2042                 faceOfLink.first = 
2043                   curLink->FirstNode()->IsLinked( curLink->LastNode()->_intPoint );
2044                 faceOfLink.second = iL;
2045               }
2046             usedFaceIDs.clear();
2047           }
2048           curFace = faceOfLink.first;
2049           curLink = freeLinks[ faceOfLink.second ];
2050           freeLinks[ faceOfLink.second ] = 0;
2051         }
2052         usedFaceIDs.insert( curFace );
2053         polygon._links.push_back( *curLink );
2054         --nbFreeLinks;
2055
2056         // find all links lying on a curFace
2057         do
2058         {
2059           // go forward from curLink
2060           curNode = curLink->LastNode();
2061           curLink = 0;
2062           for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2063             if ( freeLinks[ iL ] &&
2064                  freeLinks[ iL ]->FirstNode() == curNode &&
2065                  freeLinks[ iL ]->LastNode()->IsOnFace( curFace ))
2066             {
2067               curLink = freeLinks[ iL ];
2068               freeLinks[ iL ] = 0;
2069               polygon._links.push_back( *curLink );
2070               --nbFreeLinks;
2071             }
2072         } while ( curLink );
2073
2074         std::reverse( polygon._links.begin(), polygon._links.end() );
2075
2076         curLink = & polygon._links.back();
2077         do
2078         {
2079           // go backward from curLink
2080           curNode = curLink->FirstNode();
2081           curLink = 0;
2082           for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2083             if ( freeLinks[ iL ] &&
2084                  freeLinks[ iL ]->LastNode() == curNode &&
2085                  freeLinks[ iL ]->FirstNode()->IsOnFace( curFace ))
2086             {
2087               curLink = freeLinks[ iL ];
2088               freeLinks[ iL ] = 0;
2089               polygon._links.push_back( *curLink );
2090               --nbFreeLinks;
2091             }
2092         } while ( curLink );
2093
2094         curNode = polygon._links.back().FirstNode();
2095
2096         if ( polygon._links[0].LastNode() != curNode )
2097         {
2098           if ( nbVertexNodes > 0 )
2099           {
2100             // add links with _vIntNodes if not already used
2101             chainNodes.clear();
2102             for ( size_t iN = 0; iN < _vIntNodes.size(); ++iN )
2103               if ( !_vIntNodes[ iN ]->IsUsedInFace() &&
2104                    _vIntNodes[ iN ]->IsOnFace( curFace ))
2105               {
2106                 _vIntNodes[ iN ]->_usedInFace = &polygon;
2107                 chainNodes.push_back( _vIntNodes[ iN ] );
2108               }
2109             if ( chainNodes.size() > 1 )
2110             {
2111               sortVertexNodes( chainNodes, curNode, curFace );
2112             }
2113             for ( size_t i = 0; i < chainNodes.size(); ++i )
2114             {
2115               polygon.AddPolyLink( chainNodes[ i ], curNode );
2116               curNode = chainNodes[ i ];
2117               freeLinks.push_back( &polygon._links.back() );
2118               ++nbFreeLinks;
2119             }
2120             nbVertexNodes -= chainNodes.size();
2121           }
2122           // if ( polygon._links.size() > 1 )
2123           {
2124             polygon.AddPolyLink( polygon._links[0].LastNode(), curNode );
2125             freeLinks.push_back( &polygon._links.back() );
2126             ++nbFreeLinks;
2127           }
2128         }
2129       } // if there are intersections with EDGEs
2130
2131       if ( polygon._links.size() < 2 ||
2132            polygon._links[0].LastNode() != polygon._links.back().FirstNode() )
2133         return; // closed polygon not found -> invalid polyhedron
2134
2135       if ( polygon._links.size() == 2 )
2136       {
2137         if ( freeLinks.back() == &polygon._links.back() )
2138         {
2139           freeLinks.pop_back();
2140           --nbFreeLinks;
2141         }
2142         if ( polygon._links.front().NbFaces() > 0 )
2143           polygon._links.back().AddFace( polygon._links.front()._link->_faces[0] );
2144         if ( polygon._links.back().NbFaces() > 0 )
2145           polygon._links.front().AddFace( polygon._links.back()._link->_faces[0] );
2146
2147         if ( iPolygon == _polygons.size()-1 )
2148           _polygons.pop_back();
2149       }
2150       else // polygon._links.size() >= 2
2151       {
2152         // add polygon to its links
2153         for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
2154         {
2155           polygon._links[ iL ].AddFace( &polygon );
2156           polygon._links[ iL ].Reverse();
2157         }
2158         if ( /*hasEdgeIntersections &&*/ iPolygon == _polygons.size() - 1 )
2159         {
2160           // check that a polygon does not lie on a hexa side
2161           coplanarPolyg = 0;
2162           for ( size_t iL = 0; iL < polygon._links.size() && !coplanarPolyg; ++iL )
2163           {
2164             if ( polygon._links[ iL ].NbFaces() < 2 )
2165               continue; // it's a just added free link
2166             // look for a polygon made on a hexa side and sharing
2167             // two or more haxa links
2168             size_t iL2;
2169             coplanarPolyg = polygon._links[ iL ]._link->_faces[0];
2170             for ( iL2 = iL + 1; iL2 < polygon._links.size(); ++iL2 )
2171               if ( polygon._links[ iL2 ]._link->_faces[0] == coplanarPolyg &&
2172                    !coplanarPolyg->IsPolyLink( polygon._links[ iL  ]) &&
2173                    !coplanarPolyg->IsPolyLink( polygon._links[ iL2 ]) &&
2174                    coplanarPolyg < & _polygons[ nbQuadPolygons ])
2175                 break;
2176             if ( iL2 == polygon._links.size() )
2177               coplanarPolyg = 0;
2178           }
2179           if ( coplanarPolyg ) // coplanar polygon found
2180           {
2181             freeLinks.resize( freeLinks.size() - polygon._polyLinks.size() );
2182             nbFreeLinks -= polygon._polyLinks.size();
2183
2184             // an E_IntersectPoint used to mark nodes of coplanarPolyg
2185             // as lying on curFace while they are not at intersection with geometry
2186             ipTmp._faceIDs.resize(1);
2187             ipTmp._faceIDs[0] = curFace;
2188
2189             // fill freeLinks with links not shared by coplanarPolyg and polygon
2190             for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
2191               if ( polygon._links[ iL ]._link->_faces[1] &&
2192                    polygon._links[ iL ]._link->_faces[0] != coplanarPolyg )
2193               {
2194                 _Face* p = polygon._links[ iL ]._link->_faces[0];
2195                 for ( size_t iL2 = 0; iL2 < p->_links.size(); ++iL2 )
2196                   if ( p->_links[ iL2 ]._link == polygon._links[ iL ]._link )
2197                   {
2198                     freeLinks.push_back( & p->_links[ iL2 ] );
2199                     ++nbFreeLinks;
2200                     freeLinks.back()->RemoveFace( &polygon );
2201                     break;
2202                   }
2203               }
2204             for ( size_t iL = 0; iL < coplanarPolyg->_links.size(); ++iL )
2205               if ( coplanarPolyg->_links[ iL ]._link->_faces[1] &&
2206                    coplanarPolyg->_links[ iL ]._link->_faces[1] != &polygon )
2207               {
2208                 _Face* p = coplanarPolyg->_links[ iL ]._link->_faces[0];
2209                 if ( p == coplanarPolyg )
2210                   p = coplanarPolyg->_links[ iL ]._link->_faces[1];
2211                 for ( size_t iL2 = 0; iL2 < p->_links.size(); ++iL2 )
2212                   if ( p->_links[ iL2 ]._link == coplanarPolyg->_links[ iL ]._link )
2213                   {
2214                     // set links of coplanarPolyg in place of used freeLinks
2215                     // to re-create coplanarPolyg next
2216                     size_t iL3 = 0;
2217                     for ( ; iL3 < freeLinks.size() && freeLinks[ iL3 ]; ++iL3 );
2218                     if ( iL3 < freeLinks.size() )
2219                       freeLinks[ iL3 ] = ( & p->_links[ iL2 ] );
2220                     else
2221                       freeLinks.push_back( & p->_links[ iL2 ] );
2222                     ++nbFreeLinks;
2223                     freeLinks[ iL3 ]->RemoveFace( coplanarPolyg );
2224                     //  mark nodes of coplanarPolyg as lying on curFace
2225                     for ( int iN = 0; iN < 2; ++iN )
2226                     {
2227                       _Node* n = freeLinks[ iL3 ]->_link->_nodes[ iN ];
2228                       if ( n->_intPoint ) n->_intPoint->Add( ipTmp._faceIDs );
2229                       else                n->_intPoint = &ipTmp;
2230                     }
2231                     break;
2232                   }
2233               }
2234             // set coplanarPolyg to be re-created next
2235             for ( size_t iP = 0; iP < _polygons.size(); ++iP )
2236               if ( coplanarPolyg == & _polygons[ iP ] )
2237               {
2238                 iPolygon = iP;
2239                 _polygons[ iPolygon ]._links.clear();
2240                 _polygons[ iPolygon ]._polyLinks.clear();
2241                 break;
2242               }
2243             _polygons.pop_back();
2244             usedFaceIDs.erase( curFace );
2245             continue;
2246           } // if ( coplanarPolyg )
2247         } // if ( hasEdgeIntersections ) - search for coplanarPolyg
2248
2249         iPolygon = _polygons.size();
2250
2251       } // end of case ( polygon._links.size() > 2 )
2252     } // while ( nbFreeLinks > 0 )
2253
2254     if ( ! checkPolyhedronSize() )
2255     {
2256       return;
2257     }
2258
2259     for ( size_t i = 0; i < 8; ++i )
2260       if ( _hexNodes[ i ]._intPoint == &ipTmp )
2261         _hexNodes[ i ]._intPoint = 0;
2262
2263     // create a classic cell if possible
2264
2265     int nbPolygons = 0;
2266     for ( size_t iF = 0; iF < _polygons.size(); ++iF )
2267       nbPolygons += (_polygons[ iF ]._links.size() > 0 );
2268
2269     //const int nbNodes = _nbCornerNodes + nbIntersections;
2270     int nbNodes = 0;
2271     for ( size_t i = 0; i < 8; ++i )
2272       nbNodes += _hexNodes[ i ].IsUsedInFace();
2273     for ( size_t i = 0; i < _intNodes.size(); ++i )
2274       nbNodes += _intNodes[ i ].IsUsedInFace();
2275
2276     bool isClassicElem = false;
2277     if (      nbNodes == 8 && nbPolygons == 6 ) isClassicElem = addHexa();
2278     else if ( nbNodes == 4 && nbPolygons == 4 ) isClassicElem = addTetra();
2279     else if ( nbNodes == 6 && nbPolygons == 5 ) isClassicElem = addPenta();
2280     else if ( nbNodes == 5 && nbPolygons == 5 ) isClassicElem = addPyra ();
2281     if ( !isClassicElem )
2282     {
2283       _volumeDefs._nodes.clear();
2284       _volumeDefs._quantities.clear();
2285
2286       for ( size_t iF = 0; iF < _polygons.size(); ++iF )
2287       {
2288         const size_t nbLinks = _polygons[ iF ]._links.size();
2289         if ( nbLinks == 0 ) continue;
2290         _volumeDefs._quantities.push_back( nbLinks );
2291         for ( size_t iL = 0; iL < nbLinks; ++iL )
2292           _volumeDefs._nodes.push_back( _polygons[ iF ]._links[ iL ].FirstNode() );
2293       }
2294     }
2295   }
2296   //================================================================================
2297   /*!
2298    * \brief Create elements in the mesh
2299    */
2300   int Hexahedron::MakeElements(SMESH_MesherHelper&                      helper,
2301                                const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap)
2302   {
2303     SMESHDS_Mesh* mesh = helper.GetMeshDS();
2304
2305     size_t nbCells[3] = { _grid->_coords[0].size() - 1,
2306                           _grid->_coords[1].size() - 1,
2307                           _grid->_coords[2].size() - 1 };
2308     const size_t nbGridCells = nbCells[0] * nbCells[1] * nbCells[2];
2309     vector< Hexahedron* > allHexa( nbGridCells, 0 );
2310     int nbIntHex = 0;
2311
2312     // set intersection nodes from GridLine's to links of allHexa
2313     int i,j,k, iDirOther[3][2] = {{ 1,2 },{ 0,2 },{ 0,1 }};
2314     for ( int iDir = 0; iDir < 3; ++iDir )
2315     {
2316       int dInd[4][3] = { {0,0,0}, {0,0,0}, {0,0,0}, {0,0,0} };
2317       dInd[1][ iDirOther[iDir][0] ] = -1;
2318       dInd[2][ iDirOther[iDir][1] ] = -1;
2319       dInd[3][ iDirOther[iDir][0] ] = -1; dInd[3][ iDirOther[iDir][1] ] = -1;
2320       // loop on GridLine's parallel to iDir
2321       LineIndexer lineInd = _grid->GetLineIndexer( iDir );
2322       for ( ; lineInd.More(); ++lineInd )
2323       {
2324         GridLine& line = _grid->_lines[ iDir ][ lineInd.LineIndex() ];
2325         multiset< F_IntersectPoint >::const_iterator ip = line._intPoints.begin();
2326         for ( ; ip != line._intPoints.end(); ++ip )
2327         {
2328           // if ( !ip->_node ) continue; // intersection at a grid node
2329           lineInd.SetIndexOnLine( ip->_indexOnLine );
2330           for ( int iL = 0; iL < 4; ++iL ) // loop on 4 cells sharing a link
2331           {
2332             i = int(lineInd.I()) + dInd[iL][0];
2333             j = int(lineInd.J()) + dInd[iL][1];
2334             k = int(lineInd.K()) + dInd[iL][2];
2335             if ( i < 0 || i >= (int) nbCells[0] ||
2336                  j < 0 || j >= (int) nbCells[1] ||
2337                  k < 0 || k >= (int) nbCells[2] ) continue;
2338
2339             const size_t hexIndex = _grid->CellIndex( i,j,k );
2340             Hexahedron *& hex = allHexa[ hexIndex ];
2341             if ( !hex)
2342             {
2343               hex = new Hexahedron( *this );
2344               hex->_i = i;
2345               hex->_j = j;
2346               hex->_k = k;
2347               ++nbIntHex;
2348             }
2349             const int iLink = iL + iDir * 4;
2350             hex->_hexLinks[iLink]._fIntPoints.push_back( &(*ip) );
2351             hex->_nbFaceIntNodes += bool( ip->_node );
2352           }
2353         }
2354       }
2355     }
2356
2357     // implement geom edges into the mesh
2358     addEdges( helper, allHexa, edge2faceIDsMap );
2359
2360     // add not split hexadrons to the mesh
2361     int nbAdded = 0;
2362     vector< Hexahedron* > intHexa( nbIntHex, (Hexahedron*) NULL );
2363     for ( size_t i = 0; i < allHexa.size(); ++i )
2364     {
2365       Hexahedron * & hex = allHexa[ i ];
2366       if ( hex )
2367       {
2368         intHexa.push_back( hex );
2369         if ( hex->_nbFaceIntNodes > 0 || hex->_eIntPoints.size() > 0 )
2370           continue; // treat intersected hex later
2371         this->init( hex->_i, hex->_j, hex->_k );
2372       }
2373       else
2374       {    
2375         this->init( i );
2376       }
2377       if (( _nbCornerNodes == 8 ) &&
2378           ( _nbBndNodes < _nbCornerNodes || !isInHole() ))
2379       {
2380         // order of _hexNodes is defined by enum SMESH_Block::TShapeID
2381         SMDS_MeshElement* el =
2382           mesh->AddVolume( _hexNodes[0].Node(), _hexNodes[2].Node(),
2383                            _hexNodes[3].Node(), _hexNodes[1].Node(),
2384                            _hexNodes[4].Node(), _hexNodes[6].Node(),
2385                            _hexNodes[7].Node(), _hexNodes[5].Node() );
2386         mesh->SetMeshElementOnShape( el, helper.GetSubShapeID() );
2387         ++nbAdded;
2388         if ( hex )
2389           intHexa.pop_back();
2390       }
2391       else if ( _nbCornerNodes > 3  && !hex )
2392       {
2393         // all intersection of hex with geometry are at grid nodes
2394         hex = new Hexahedron( *this );
2395         hex->_i = _i;
2396         hex->_j = _j;
2397         hex->_k = _k;
2398         intHexa.push_back( hex );
2399       }
2400     }
2401
2402     // add elements resulted from hexadron intersection
2403 #ifdef WITH_TBB
2404     tbb::parallel_for ( tbb::blocked_range<size_t>( 0, intHexa.size() ),
2405                         ParallelHexahedron( intHexa ),
2406                         tbb::simple_partitioner()); // ComputeElements() is called here
2407     for ( size_t i = 0; i < intHexa.size(); ++i )
2408       if ( Hexahedron * hex = intHexa[ i ] )
2409         nbAdded += hex->addElements( helper );
2410 #else
2411     for ( size_t i = 0; i < intHexa.size(); ++i )
2412       if ( Hexahedron * hex = intHexa[ i ] )
2413       {
2414         hex->ComputeElements();
2415         nbAdded += hex->addElements( helper );
2416       }
2417 #endif
2418
2419     for ( size_t i = 0; i < allHexa.size(); ++i )
2420       if ( allHexa[ i ] )
2421         delete allHexa[ i ];
2422
2423     return nbAdded;
2424   }
2425
2426   //================================================================================
2427   /*!
2428    * \brief Implements geom edges into the mesh
2429    */
2430   void Hexahedron::addEdges(SMESH_MesherHelper&                      helper,
2431                             vector< Hexahedron* >&                   hexes,
2432                             const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap)
2433   {
2434     if ( edge2faceIDsMap.empty() ) return;
2435
2436     // Prepare planes for intersecting with EDGEs
2437     GridPlanes pln[3];
2438     {
2439       for ( int iDirZ = 0; iDirZ < 3; ++iDirZ ) // iDirZ gives normal direction to planes
2440       {
2441         GridPlanes& planes = pln[ iDirZ ];
2442         int iDirX = ( iDirZ + 1 ) % 3;
2443         int iDirY = ( iDirZ + 2 ) % 3;
2444         planes._zNorm  = ( _grid->_axes[ iDirX ] ^ _grid->_axes[ iDirY ] ).Normalized();
2445         planes._zProjs.resize ( _grid->_coords[ iDirZ ].size() );
2446         planes._zProjs [0] = 0;
2447         const double       zFactor = _grid->_axes[ iDirZ ] * planes._zNorm;
2448         const vector< double > & u = _grid->_coords[ iDirZ ];
2449         for ( size_t i = 1; i < planes._zProjs.size(); ++i )
2450         {
2451           planes._zProjs [i] = zFactor * ( u[i] - u[0] );
2452         }
2453       }
2454     }
2455     const double deflection = _grid->_minCellSize / 20.;
2456     const double tol        = _grid->_tol;
2457     E_IntersectPoint ip;
2458
2459     // Intersect EDGEs with the planes
2460     map< TGeomID, vector< TGeomID > >::const_iterator e2fIt = edge2faceIDsMap.begin();
2461     for ( ; e2fIt != edge2faceIDsMap.end(); ++e2fIt )
2462     {
2463       const TGeomID  edgeID = e2fIt->first;
2464       const TopoDS_Edge & E = TopoDS::Edge( _grid->_shapes( edgeID ));
2465       BRepAdaptor_Curve curve( E );
2466       TopoDS_Vertex v1 = helper.IthVertex( 0, E, false ); 
2467       TopoDS_Vertex v2 = helper.IthVertex( 1, E, false ); 
2468
2469       ip._faceIDs = e2fIt->second;
2470       ip._shapeID = edgeID;
2471
2472       // discretize the EDGE
2473       GCPnts_UniformDeflection discret( curve, deflection, true );
2474       if ( !discret.IsDone() || discret.NbPoints() < 2 )
2475         continue;
2476
2477       // perform intersection
2478       for ( int iDirZ = 0; iDirZ < 3; ++iDirZ )
2479       {
2480         GridPlanes& planes = pln[ iDirZ ];
2481         int      iDirX = ( iDirZ + 1 ) % 3;
2482         int      iDirY = ( iDirZ + 2 ) % 3;
2483         double    xLen = _grid->_coords[ iDirX ].back() - _grid->_coords[ iDirX ][0];
2484         double    yLen = _grid->_coords[ iDirY ].back() - _grid->_coords[ iDirY ][0];
2485         double    zLen = _grid->_coords[ iDirZ ].back() - _grid->_coords[ iDirZ ][0];
2486         int dIJK[3], d000[3] = { 0,0,0 };
2487         double o[3] = { _grid->_coords[0][0],
2488                         _grid->_coords[1][0],
2489                         _grid->_coords[2][0] };
2490
2491         // locate the 1st point of a segment within the grid
2492         gp_XYZ p1     = discret.Value( 1 ).XYZ();
2493         double u1     = discret.Parameter( 1 );
2494         double zProj1 = planes._zNorm * ( p1 - _grid->_origin );
2495
2496         _grid->ComputeUVW( p1, ip._uvw );
2497         int iX1 = int(( ip._uvw[iDirX] - o[iDirX]) / xLen * (_grid->_coords[ iDirX ].size() - 1));
2498         int iY1 = int(( ip._uvw[iDirY] - o[iDirY]) / yLen * (_grid->_coords[ iDirY ].size() - 1));
2499         int iZ1 = int(( ip._uvw[iDirZ] - o[iDirZ]) / zLen * (_grid->_coords[ iDirZ ].size() - 1));
2500         locateValue( iX1, ip._uvw[iDirX], _grid->_coords[ iDirX ], dIJK[ iDirX ], tol );
2501         locateValue( iY1, ip._uvw[iDirY], _grid->_coords[ iDirY ], dIJK[ iDirY ], tol );
2502         locateValue( iZ1, ip._uvw[iDirZ], _grid->_coords[ iDirZ ], dIJK[ iDirZ ], tol );
2503
2504         int ijk[3]; // grid index where a segment intersect a plane
2505         ijk[ iDirX ] = iX1;
2506         ijk[ iDirY ] = iY1;
2507         ijk[ iDirZ ] = iZ1;
2508
2509         // add the 1st vertex point to a hexahedron
2510         if ( iDirZ == 0 )
2511         {
2512           ip._point   = p1;
2513           ip._shapeID = _grid->_shapes.Add( v1 );
2514           _grid->_edgeIntP.push_back( ip );
2515           if ( !addIntersection( _grid->_edgeIntP.back(), hexes, ijk, d000 ))
2516             _grid->_edgeIntP.pop_back();
2517           ip._shapeID = edgeID;
2518         }
2519         for ( int iP = 2; iP <= discret.NbPoints(); ++iP )
2520         {
2521           // locate the 2nd point of a segment within the grid
2522           gp_XYZ p2     = discret.Value( iP ).XYZ();
2523           double u2     = discret.Parameter( iP );
2524           double zProj2 = planes._zNorm * ( p2 - _grid->_origin );
2525           int    iZ2    = iZ1;
2526           if ( Abs( zProj2 - zProj1 ) > std::numeric_limits<double>::min() )
2527           {
2528             locateValue( iZ2, zProj2, planes._zProjs, dIJK[ iDirZ ], tol );
2529
2530             // treat intersections with planes between 2 end points of a segment
2531             int dZ = ( iZ1 <= iZ2 ) ? +1 : -1;
2532             int iZ = iZ1 + ( iZ1 < iZ2 );
2533             for ( int i = 0, nb = Abs( iZ1 - iZ2 ); i < nb; ++i, iZ += dZ )
2534             {
2535               ip._point = findIntPoint( u1, zProj1, u2, zProj2,
2536                                         planes._zProjs[ iZ ],
2537                                         curve, planes._zNorm, _grid->_origin );
2538               _grid->ComputeUVW( ip._point.XYZ(), ip._uvw );
2539               locateValue( ijk[iDirX], ip._uvw[iDirX], _grid->_coords[iDirX], dIJK[iDirX], tol );
2540               locateValue( ijk[iDirY], ip._uvw[iDirY], _grid->_coords[iDirY], dIJK[iDirY], tol );
2541               ijk[ iDirZ ] = iZ;
2542
2543               // add ip to hex "above" the plane
2544               _grid->_edgeIntP.push_back( ip );
2545               dIJK[ iDirZ ] = 0;
2546               bool added = addIntersection(_grid->_edgeIntP.back(), hexes, ijk, dIJK);
2547
2548               // add ip to hex "below" the plane
2549               ijk[ iDirZ ] = iZ-1;
2550               if ( !addIntersection( _grid->_edgeIntP.back(), hexes, ijk, dIJK ) &&
2551                    !added)
2552                 _grid->_edgeIntP.pop_back();
2553             }
2554           }
2555           iZ1    = iZ2;
2556           p1     = p2;
2557           u1     = u2;
2558           zProj1 = zProj2;
2559         }
2560         // add the 2nd vertex point to a hexahedron
2561         if ( iDirZ == 0 )
2562         {
2563           ip._shapeID = _grid->_shapes.Add( v2 );
2564           ip._point = p1;
2565           _grid->ComputeUVW( p1, ip._uvw );
2566           locateValue( ijk[iDirX], ip._uvw[iDirX], _grid->_coords[iDirX], dIJK[iDirX], tol );
2567           locateValue( ijk[iDirY], ip._uvw[iDirY], _grid->_coords[iDirY], dIJK[iDirY], tol );
2568           ijk[ iDirZ ] = iZ1;
2569           _grid->_edgeIntP.push_back( ip );
2570           if ( !addIntersection( _grid->_edgeIntP.back(), hexes, ijk, d000 ))
2571             _grid->_edgeIntP.pop_back();
2572           ip._shapeID = edgeID;
2573         }
2574       } // loop on 3 grid directions
2575     } // loop on EDGEs
2576
2577   }
2578
2579   //================================================================================
2580   /*!
2581    * \brief Finds intersection of a curve with a plane
2582    *  \param [in] u1 - parameter of one curve point
2583    *  \param [in] proj1 - projection of the curve point to the plane normal
2584    *  \param [in] u2 - parameter of another curve point
2585    *  \param [in] proj2 - projection of the other curve point to the plane normal
2586    *  \param [in] proj - projection of a point where the curve intersects the plane
2587    *  \param [in] curve - the curve
2588    *  \param [in] axis - the plane normal
2589    *  \param [in] origin - the plane origin
2590    *  \return gp_Pnt - the found intersection point
2591    */
2592   gp_Pnt Hexahedron::findIntPoint( double u1, double proj1,
2593                                    double u2, double proj2,
2594                                    double proj,
2595                                    BRepAdaptor_Curve& curve,
2596                                    const gp_XYZ& axis,
2597                                    const gp_XYZ& origin)
2598   {
2599     double r = (( proj - proj1 ) / ( proj2 - proj1 ));
2600     double u = u1 * ( 1 - r ) + u2 * r;
2601     gp_Pnt p = curve.Value( u );
2602     double newProj =  axis * ( p.XYZ() - origin );
2603     if ( Abs( proj - newProj ) > _grid->_tol / 10. )
2604     {
2605       if ( r > 0.5 )
2606         return findIntPoint( u2, proj2, u, newProj, proj, curve, axis, origin );
2607       else
2608         return findIntPoint( u1, proj2, u, newProj, proj, curve, axis, origin );
2609     }
2610     return p;
2611   }
2612
2613   //================================================================================
2614   /*!
2615    * \brief Returns indices of a hexahedron sub-entities holding a point
2616    *  \param [in] ip - intersection point
2617    *  \param [out] facets - 0-3 facets holding a point
2618    *  \param [out] sub - index of a vertex or an edge holding a point
2619    *  \return int - number of facets holding a point
2620    */
2621   int Hexahedron::getEntity( const E_IntersectPoint* ip, int* facets, int& sub )
2622   {
2623     enum { X = 1, Y = 2, Z = 4 }; // == 001, 010, 100
2624     int nbFacets = 0;
2625     int vertex = 0, edgeMask = 0;
2626
2627     if ( Abs( _grid->_coords[0][ _i   ] - ip->_uvw[0] ) < _grid->_tol ) {
2628       facets[ nbFacets++ ] = SMESH_Block::ID_F0yz;
2629       edgeMask |= X;
2630     }
2631     else if ( Abs( _grid->_coords[0][ _i+1 ] - ip->_uvw[0] ) < _grid->_tol ) {
2632       facets[ nbFacets++ ] = SMESH_Block::ID_F1yz;
2633       vertex   |= X;
2634       edgeMask |= X;
2635     }
2636     if ( Abs( _grid->_coords[1][ _j   ] - ip->_uvw[1] ) < _grid->_tol ) {
2637       facets[ nbFacets++ ] = SMESH_Block::ID_Fx0z;
2638       edgeMask |= Y;
2639     }
2640     else if ( Abs( _grid->_coords[1][ _j+1 ] - ip->_uvw[1] ) < _grid->_tol ) {
2641       facets[ nbFacets++ ] = SMESH_Block::ID_Fx1z;
2642       vertex   |= Y;
2643       edgeMask |= Y;
2644     }
2645     if ( Abs( _grid->_coords[2][ _k   ] - ip->_uvw[2] ) < _grid->_tol ) {
2646       facets[ nbFacets++ ] = SMESH_Block::ID_Fxy0;
2647       edgeMask |= Z;
2648     }
2649     else if ( Abs( _grid->_coords[2][ _k+1 ] - ip->_uvw[2] ) < _grid->_tol ) {
2650       facets[ nbFacets++ ] = SMESH_Block::ID_Fxy1;
2651       vertex   |= Z;
2652       edgeMask |= Z;
2653     }
2654
2655     switch ( nbFacets )
2656     {
2657     case 0: sub = 0;         break;
2658     case 1: sub = facets[0]; break;
2659     case 2: {
2660       const int edge [3][8] = {
2661         { SMESH_Block::ID_E00z, SMESH_Block::ID_E10z,
2662           SMESH_Block::ID_E01z, SMESH_Block::ID_E11z },
2663         { SMESH_Block::ID_E0y0, SMESH_Block::ID_E1y0, 0, 0,
2664           SMESH_Block::ID_E0y1, SMESH_Block::ID_E1y1 },
2665         { SMESH_Block::ID_Ex00, 0, SMESH_Block::ID_Ex10, 0,
2666           SMESH_Block::ID_Ex01, 0, SMESH_Block::ID_Ex11 }
2667       };
2668       switch ( edgeMask ) {
2669       case X | Y: sub = edge[ 0 ][ vertex ]; break;
2670       case X | Z: sub = edge[ 1 ][ vertex ]; break;
2671       default:    sub = edge[ 2 ][ vertex ];
2672       }
2673       break;
2674     }
2675     //case 3:
2676     default:
2677       sub = vertex + SMESH_Block::ID_FirstV;
2678     }
2679
2680     return nbFacets;
2681   }
2682   //================================================================================
2683   /*!
2684    * \brief Adds intersection with an EDGE
2685    */
2686   bool Hexahedron::addIntersection( const E_IntersectPoint& ip,
2687                                     vector< Hexahedron* >&  hexes,
2688                                     int ijk[], int dIJK[] )
2689   {
2690     bool added = false;
2691
2692     size_t hexIndex[4] = {
2693       _grid->CellIndex( ijk[0], ijk[1], ijk[2] ),
2694       dIJK[0] ? _grid->CellIndex( ijk[0]+dIJK[0], ijk[1], ijk[2] ) : -1,
2695       dIJK[1] ? _grid->CellIndex( ijk[0], ijk[1]+dIJK[1], ijk[2] ) : -1,
2696       dIJK[2] ? _grid->CellIndex( ijk[0], ijk[1], ijk[2]+dIJK[2] ) : -1
2697     };
2698     for ( int i = 0; i < 4; ++i )
2699     {
2700       if ( /*0 <= hexIndex[i] &&*/ hexIndex[i] < hexes.size() && hexes[ hexIndex[i] ] )
2701       {
2702         Hexahedron* h = hexes[ hexIndex[i] ];
2703         // check if ip is really inside the hex
2704 #ifdef _DEBUG_
2705         if ( h->isOutParam( ip._uvw ))
2706           throw SALOME_Exception("ip outside a hex");
2707 #endif
2708         h->_eIntPoints.push_back( & ip );
2709         added = true;
2710       }
2711     }
2712     return added;
2713   }
2714   //================================================================================
2715   /*!
2716    * \brief Finds nodes at a path from one node to another via intersections with EDGEs
2717    */
2718   bool Hexahedron::findChain( _Node*          n1,
2719                               _Node*          n2,
2720                               _Face&          quad,
2721                               vector<_Node*>& chn )
2722   {
2723     chn.clear();
2724     chn.push_back( n1 );
2725     for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
2726       if ( !quad._eIntNodes[ iP ]->IsUsedInFace( &quad ) &&
2727            n1->IsLinked( quad._eIntNodes[ iP ]->_intPoint ) &&
2728            n2->IsLinked( quad._eIntNodes[ iP ]->_intPoint ))
2729       {
2730         chn.push_back( quad._eIntNodes[ iP ]);
2731         chn.push_back( n2 );
2732         quad._eIntNodes[ iP ]->_usedInFace = &quad;
2733         return true;
2734       }
2735     bool found;
2736     do
2737     {
2738       found = false;
2739       for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
2740         if ( !quad._eIntNodes[ iP ]->IsUsedInFace( &quad ) &&
2741              chn.back()->IsLinked( quad._eIntNodes[ iP ]->_intPoint ))
2742         {
2743           chn.push_back( quad._eIntNodes[ iP ]);
2744           found = ( quad._eIntNodes[ iP ]->_usedInFace = &quad );
2745           break;
2746         }
2747     } while ( found && ! chn.back()->IsLinked( n2->_intPoint ) );
2748
2749     if ( chn.back() != n2 && chn.back()->IsLinked( n2->_intPoint ))
2750       chn.push_back( n2 );
2751
2752     return chn.size() > 1;
2753   }
2754   //================================================================================
2755   /*!
2756    * \brief Try to heal a polygon whose ends are not connected
2757    */
2758   bool Hexahedron::closePolygon( _Face* polygon, vector<_Node*>& chainNodes ) const
2759   {
2760     int i = -1, nbLinks = polygon->_links.size();
2761     if ( nbLinks < 3 )
2762       return false;
2763     vector< _OrientedLink > newLinks;
2764     // find a node lying on the same FACE as the last one
2765     _Node*   node = polygon->_links.back().LastNode();
2766     int avoidFace = node->IsLinked( polygon->_links.back().FirstNode()->_intPoint );
2767     for ( i = nbLinks - 2; i >= 0; --i )
2768       if ( node->IsLinked( polygon->_links[i].FirstNode()->_intPoint, avoidFace ))
2769         break;
2770     if ( i >= 0 )
2771     {
2772       for ( ; i < nbLinks; ++i )
2773         newLinks.push_back( polygon->_links[i] );
2774     }
2775     else
2776     {
2777       // find a node lying on the same FACE as the first one
2778       node      = polygon->_links[0].FirstNode();
2779       avoidFace = node->IsLinked( polygon->_links[0].LastNode()->_intPoint );
2780       for ( i = 1; i < nbLinks; ++i )
2781         if ( node->IsLinked( polygon->_links[i].LastNode()->_intPoint, avoidFace ))
2782           break;
2783       if ( i < nbLinks )
2784         for ( nbLinks = i + 1, i = 0; i < nbLinks; ++i )
2785           newLinks.push_back( polygon->_links[i] );
2786     }
2787     if ( newLinks.size() > 1 )
2788     {
2789       polygon->_links.swap( newLinks );
2790       chainNodes.clear();
2791       chainNodes.push_back( polygon->_links.back().LastNode() );
2792       chainNodes.push_back( polygon->_links[0].FirstNode() );
2793       return true;
2794     }
2795     return false;
2796   }
2797   //================================================================================
2798   /*!
2799    * \brief Finds nodes on the same EDGE as the first node of avoidSplit.
2800    *
2801    * This function is for a case where an EDGE lies on a quad which lies on a FACE
2802    * so that a part of quad in ON and another part in IN
2803    */
2804   bool Hexahedron::findChainOnEdge( const vector< _OrientedLink >& splits,
2805                                     const _OrientedLink&           prevSplit,
2806                                     const _OrientedLink&           avoidSplit,
2807                                     size_t &                       iS,
2808                                     _Face&                         quad,
2809                                     vector<_Node*>&                chn )
2810   {
2811     if ( !isImplementEdges() )
2812       return false;
2813
2814     _Node* pn1 = prevSplit.FirstNode();
2815     _Node* pn2 = prevSplit.LastNode();
2816     int avoidFace = pn1->IsLinked( pn2->_intPoint ); // FACE under the quad
2817     if ( avoidFace < 1 && pn1->_intPoint )
2818       return false;
2819
2820     _Node* n, *stopNode = avoidSplit.LastNode();
2821
2822     chn.clear();
2823     if ( !quad._eIntNodes.empty() )
2824     {
2825       chn.push_back( pn2 );
2826       bool found;
2827       do
2828       {
2829         found = false;
2830         for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
2831           if (( !quad._eIntNodes[ iP ]->IsUsedInFace( &quad )) &&
2832               ( chn.back()->IsLinked( quad._eIntNodes[ iP ]->_intPoint, avoidFace )) &&
2833               ( !avoidFace || quad._eIntNodes[ iP ]->IsOnFace( avoidFace )))
2834           {
2835             chn.push_back( quad._eIntNodes[ iP ]);
2836             found = ( quad._eIntNodes[ iP ]->_usedInFace = &quad );
2837             break;
2838           }
2839       } while ( found );
2840       pn2 = chn.back();
2841     }
2842
2843     int i;
2844     for ( i = splits.size()-1; i >= 0; --i )
2845     {
2846       if ( !splits[i] )
2847         continue;
2848
2849       n = splits[i].LastNode();
2850       if ( n == stopNode )
2851         break;
2852       if (( n != pn1 ) &&
2853           ( n->IsLinked( pn2->_intPoint, avoidFace )) &&
2854           ( !avoidFace || n->IsOnFace( avoidFace )))
2855         break;
2856
2857       n = splits[i].FirstNode();
2858       if ( n == stopNode )
2859         break;
2860       if (( n->IsLinked( pn2->_intPoint, avoidFace )) &&
2861           ( !avoidFace || n->IsOnFace( avoidFace )))
2862         break;
2863       n = 0;
2864     }
2865     if ( n && n != stopNode)
2866     {
2867       if ( chn.empty() )
2868         chn.push_back( pn2 );
2869       chn.push_back( n );
2870       iS = i-1;
2871       return true;
2872     }
2873     return false;
2874   }
2875   //================================================================================
2876   /*!
2877    * \brief Checks transition at the ginen intersection node of a link
2878    */
2879   bool Hexahedron::isOutPoint( _Link& link, int iP, SMESH_MesherHelper& helper ) const
2880   {
2881     bool isOut = false;
2882
2883     const bool moreIntPoints = ( iP+1 < (int) link._fIntPoints.size() );
2884
2885     // get 2 _Node's
2886     _Node* n1 = link._fIntNodes[ iP ];
2887     if ( !n1->Node() )
2888       n1 = link._nodes[0];
2889     _Node* n2 = moreIntPoints ? link._fIntNodes[ iP+1 ] : 0;
2890     if ( !n2 || !n2->Node() )
2891       n2 = link._nodes[1];
2892     if ( !n2->Node() )
2893       return true;
2894
2895     // get all FACEs under n1 and n2
2896     set< TGeomID > faceIDs;
2897     if ( moreIntPoints ) faceIDs.insert( link._fIntPoints[iP+1]->_faceIDs.begin(),
2898                                          link._fIntPoints[iP+1]->_faceIDs.end() );
2899     if ( n2->_intPoint ) faceIDs.insert( n2->_intPoint->_faceIDs.begin(),
2900                                          n2->_intPoint->_faceIDs.end() );
2901     if ( faceIDs.empty() )
2902       return false; // n2 is inside
2903     if ( n1->_intPoint ) faceIDs.insert( n1->_intPoint->_faceIDs.begin(),
2904                                          n1->_intPoint->_faceIDs.end() );
2905     faceIDs.insert( link._fIntPoints[iP]->_faceIDs.begin(),
2906                     link._fIntPoints[iP]->_faceIDs.end() );
2907
2908     // get a point between 2 nodes
2909     gp_Pnt p1      = n1->Point();
2910     gp_Pnt p2      = n2->Point();
2911     gp_Pnt pOnLink = 0.8 * p1.XYZ() + 0.2 * p2.XYZ();
2912
2913     TopLoc_Location loc;
2914
2915     set< TGeomID >::iterator faceID = faceIDs.begin();
2916     for ( ; faceID != faceIDs.end(); ++faceID )
2917     {
2918       // project pOnLink on a FACE
2919       if ( *faceID < 1 ) continue;
2920       const TopoDS_Face& face = TopoDS::Face( _grid->_shapes( *faceID ));
2921       GeomAPI_ProjectPointOnSurf& proj =
2922         helper.GetProjector( face, loc, 0.1*_grid->_tol );
2923       gp_Pnt testPnt = pOnLink.Transformed( loc.Transformation().Inverted() );
2924       proj.Perform( testPnt );
2925       if ( proj.IsDone() && proj.NbPoints() > 0 )       
2926       {
2927         Standard_Real u,v;
2928         proj.LowerDistanceParameters( u,v );
2929
2930         if ( proj.LowerDistance() <= 0.1 * _grid->_tol )
2931         {
2932           isOut = false;
2933         }
2934         else
2935         {
2936           // find isOut by normals
2937           gp_Dir normal;
2938           if ( GeomLib::NormEstim( BRep_Tool::Surface( face, loc ),
2939                                    gp_Pnt2d( u,v ),
2940                                    0.1*_grid->_tol,
2941                                    normal ) < 3 )
2942           {
2943             if ( face.Orientation() == TopAbs_REVERSED )
2944               normal.Reverse();
2945             gp_Vec v( proj.NearestPoint(), testPnt );
2946             isOut = ( v * normal > 0 );
2947           }
2948         }
2949         if ( !isOut )
2950         {
2951           // classify a projection
2952           if ( !n1->IsOnFace( *faceID ) || !n2->IsOnFace( *faceID ))
2953           {
2954             BRepTopAdaptor_FClass2d cls( face, Precision::Confusion() );
2955             TopAbs_State state = cls.Perform( gp_Pnt2d( u,v ));
2956             if ( state == TopAbs_OUT )
2957             {
2958               isOut = true;
2959               continue;
2960             }
2961           }
2962           return false;
2963         }
2964       }
2965     }
2966     return isOut;
2967   }
2968   //================================================================================
2969   /*!
2970    * \brief Sort nodes on a FACE
2971    */
2972   void Hexahedron::sortVertexNodes(vector<_Node*>& nodes, _Node* curNode, TGeomID faceID)
2973   {
2974     if ( nodes.size() > 20 ) return;
2975
2976     // get shapes under nodes
2977     TGeomID nShapeIds[20], *nShapeIdsEnd = &nShapeIds[0] + nodes.size();
2978     for ( size_t i = 0; i < nodes.size(); ++i )
2979       if ( !( nShapeIds[i] = nodes[i]->ShapeID() ))
2980         return;
2981
2982     // get shapes of the FACE
2983     const TopoDS_Face&  face = TopoDS::Face( _grid->_shapes( faceID ));
2984     list< TopoDS_Edge > edges;
2985     list< int >         nbEdges;
2986     int nbW = SMESH_Block::GetOrderedEdges (face, edges, nbEdges);
2987     if ( nbW > 1 ) {
2988       // select a WIRE - remove EDGEs of irrelevant WIREs from edges
2989       list< TopoDS_Edge >::iterator e = edges.begin(), eEnd = e;
2990       list< int >::iterator nE = nbEdges.begin();
2991       for ( ; nbW > 0; ++nE, --nbW )
2992       {
2993         std::advance( eEnd, *nE );
2994         for ( ; e != eEnd; ++e )
2995           for ( int i = 0; i < 2; ++i )
2996           {
2997             TGeomID id = i==0 ?
2998               _grid->_shapes.FindIndex( *e ) :
2999               _grid->_shapes.FindIndex( SMESH_MesherHelper::IthVertex( 0, *e ));
3000             if (( id > 0 ) &&
3001                 ( std::find( &nShapeIds[0], nShapeIdsEnd, id ) != nShapeIdsEnd ))
3002             {
3003               edges.erase( eEnd, edges.end() ); // remove rest wires
3004               e = eEnd = edges.end();
3005               --e;
3006               nbW = 0;
3007               break;
3008             }
3009           }
3010         if ( nbW > 0 )
3011           edges.erase( edges.begin(), eEnd ); // remove a current irrelevant wire
3012       }
3013     }
3014     // rotate edges to have the first one at least partially out of the hexa
3015     list< TopoDS_Edge >::iterator e = edges.begin(), eMidOut = edges.end();
3016     for ( ; e != edges.end(); ++e )
3017     {
3018       if ( !_grid->_shapes.FindIndex( *e ))
3019         continue;
3020       bool isOut = false;
3021       gp_Pnt p;
3022       double uvw[3], f,l;
3023       for ( int i = 0; i < 2 && !isOut; ++i )
3024       {
3025         if ( i == 0 )
3026         {
3027           TopoDS_Vertex v = SMESH_MesherHelper::IthVertex( 0, *e );
3028           p = BRep_Tool::Pnt( v );
3029         }
3030         else if ( eMidOut == edges.end() )
3031         {
3032           TopLoc_Location loc;
3033           Handle(Geom_Curve) c = BRep_Tool::Curve( *e, loc, f, l);
3034           if ( c.IsNull() ) break;
3035           p = c->Value( 0.5 * ( f + l )).Transformed( loc );
3036         }
3037         else
3038         {
3039           continue;
3040         }
3041
3042         _grid->ComputeUVW( p.XYZ(), uvw );
3043         if ( isOutParam( uvw ))
3044         {
3045           if ( i == 0 )
3046             isOut = true;
3047           else
3048             eMidOut = e;
3049         }
3050       }
3051       if ( isOut )
3052         break;
3053     }
3054     if ( e != edges.end() )
3055       edges.splice( edges.end(), edges, edges.begin(), e );
3056     else if ( eMidOut != edges.end() )
3057       edges.splice( edges.end(), edges, edges.begin(), eMidOut );
3058
3059     // sort nodes according to the order of edges
3060     _Node*  orderNodes   [20];
3061     //TGeomID orderShapeIDs[20];
3062     size_t nbN = 0;
3063     TGeomID id, *pID = 0;
3064     for ( e = edges.begin(); e != edges.end(); ++e )
3065     {
3066       if (( id = _grid->_shapes.FindIndex( SMESH_MesherHelper::IthVertex( 0, *e ))) &&
3067           (( pID = std::find( &nShapeIds[0], nShapeIdsEnd, id )) != nShapeIdsEnd ))
3068       {
3069         //orderShapeIDs[ nbN ] = id;
3070         orderNodes   [ nbN++ ] = nodes[ pID - &nShapeIds[0] ];
3071         *pID = -1;
3072       }
3073       if (( id = _grid->_shapes.FindIndex( *e )) &&
3074           (( pID = std::find( &nShapeIds[0], nShapeIdsEnd, id )) != nShapeIdsEnd ))
3075       {
3076         //orderShapeIDs[ nbN ] = id;
3077         orderNodes   [ nbN++ ] = nodes[ pID - &nShapeIds[0] ];
3078         *pID = -1;
3079       }
3080     }
3081     if ( nbN != nodes.size() )
3082       return;
3083
3084     bool reverse = ( orderNodes[0    ]->Point().SquareDistance( curNode->Point() ) >
3085                      orderNodes[nbN-1]->Point().SquareDistance( curNode->Point() ));
3086
3087     for ( size_t i = 0; i < nodes.size(); ++i )
3088       nodes[ i ] = orderNodes[ reverse ? nbN-1-i : i ];
3089   }
3090
3091   //================================================================================
3092   /*!
3093    * \brief Adds computed elements to the mesh
3094    */
3095   int Hexahedron::addElements(SMESH_MesherHelper& helper)
3096   {
3097     int nbAdded = 0;
3098     // add elements resulted from hexahedron intersection
3099     //for ( size_t i = 0; i < _volumeDefs.size(); ++i )
3100     {
3101       vector< const SMDS_MeshNode* > nodes( _volumeDefs._nodes.size() );
3102       for ( size_t iN = 0; iN < nodes.size(); ++iN )
3103         if ( !( nodes[iN] = _volumeDefs._nodes[iN]->Node() ))
3104         {
3105           if ( const E_IntersectPoint* eip = _volumeDefs._nodes[iN]->EdgeIntPnt() )
3106             nodes[iN] = _volumeDefs._nodes[iN]->_intPoint->_node =
3107               helper.AddNode( eip->_point.X(),
3108                               eip->_point.Y(),
3109                               eip->_point.Z() );
3110           else
3111             throw SALOME_Exception("Bug: no node at intersection point");
3112         }
3113
3114       if ( !_volumeDefs._quantities.empty() )
3115       {
3116         helper.AddPolyhedralVolume( nodes, _volumeDefs._quantities );
3117       }
3118       else
3119       {
3120         switch ( nodes.size() )
3121         {
3122         case 8: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],
3123                                   nodes[4],nodes[5],nodes[6],nodes[7] );
3124           break;
3125         case 4: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3] );
3126           break;
3127         case 6: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3], nodes[4],nodes[5] );
3128           break;
3129         case 5:
3130           helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],nodes[4] );
3131           break;
3132         }
3133       }
3134       nbAdded += int ( _volumeDefs._nodes.size() > 0 );
3135     }
3136
3137     return nbAdded;
3138   }
3139   //================================================================================
3140   /*!
3141    * \brief Return true if the element is in a hole
3142    */
3143   bool Hexahedron::isInHole() const
3144   {
3145     if ( !_vIntNodes.empty() )
3146       return false;
3147
3148     const size_t ijk[3] = { _i, _j, _k };
3149     F_IntersectPoint curIntPnt;
3150
3151     // consider a cell to be in a hole if all links in any direction
3152     // comes OUT of geometry
3153     for ( int iDir = 0; iDir < 3; ++iDir )
3154     {
3155       const vector<double>& coords = _grid->_coords[ iDir ];
3156       LineIndexer               li = _grid->GetLineIndexer( iDir );
3157       li.SetIJK( _i,_j,_k );
3158       size_t lineIndex[4] = { li.LineIndex  (),
3159                               li.LineIndex10(),
3160                               li.LineIndex01(),
3161                               li.LineIndex11() };
3162       bool allLinksOut = true, hasLinks = false;
3163       for ( int iL = 0; iL < 4 && allLinksOut; ++iL ) // loop on 4 links parallel to iDir
3164       {
3165         const _Link& link = _hexLinks[ iL + 4*iDir ];
3166         // check transition of the first node of a link
3167         const F_IntersectPoint* firstIntPnt = 0;
3168         if ( link._nodes[0]->Node() ) // 1st node is a hexa corner
3169         {
3170           curIntPnt._paramOnLine = coords[ ijk[ iDir ]] - coords[0] + _grid->_tol;
3171           const GridLine& line = _grid->_lines[ iDir ][ lineIndex[ iL ]];
3172           multiset< F_IntersectPoint >::const_iterator ip =
3173             line._intPoints.upper_bound( curIntPnt );
3174           --ip;
3175           firstIntPnt = &(*ip);
3176         }
3177         else if ( !link._fIntPoints.empty() )
3178         {
3179           firstIntPnt = link._fIntPoints[0];
3180         }
3181
3182         if ( firstIntPnt )
3183         {
3184           hasLinks = true;
3185           allLinksOut = ( firstIntPnt->_transition == Trans_OUT );
3186         }
3187       }
3188       if ( hasLinks && allLinksOut )
3189         return true;
3190     }
3191     return false;
3192   }
3193
3194   //================================================================================
3195   /*!
3196    * \brief Return true if a polyhedron passes _sizeThreshold criterion
3197    */
3198   bool Hexahedron::checkPolyhedronSize() const
3199   {
3200     double volume = 0;
3201     for ( size_t iP = 0; iP < _polygons.size(); ++iP )
3202     {
3203       const _Face& polygon = _polygons[iP];
3204       if ( polygon._links.empty() )
3205         continue;
3206       gp_XYZ area (0,0,0);
3207       gp_XYZ p1 = polygon._links[ 0 ].FirstNode()->Point().XYZ();
3208       for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
3209       {
3210         gp_XYZ p2 = polygon._links[ iL ].LastNode()->Point().XYZ();
3211         area += p1 ^ p2;
3212         p1 = p2;
3213       }
3214       volume += p1 * area;
3215     }
3216     volume /= 6;
3217
3218     double initVolume = _sideLength[0] * _sideLength[1] * _sideLength[2];
3219
3220     return volume > initVolume / _sizeThreshold;
3221   }
3222   //================================================================================
3223   /*!
3224    * \brief Tries to create a hexahedron
3225    */
3226   bool Hexahedron::addHexa()
3227   {
3228     int nbQuad = 0, iQuad = -1;
3229     for ( size_t i = 0; i < _polygons.size(); ++i )
3230     {
3231       if ( _polygons[i]._links.empty() )
3232         continue;
3233       if ( _polygons[i]._links.size() != 4 )
3234         return false;
3235       ++nbQuad;
3236       if ( iQuad < 0 )
3237         iQuad = i;
3238     }
3239     if ( nbQuad != 6 )
3240       return false;
3241
3242     _Node* nodes[8];
3243     int nbN = 0;
3244     for ( int iL = 0; iL < 4; ++iL )
3245     {
3246       // a base node
3247       nodes[iL] = _polygons[iQuad]._links[iL].FirstNode();
3248       ++nbN;
3249
3250       // find a top node above the base node
3251       _Link* link = _polygons[iQuad]._links[iL]._link;
3252       if ( !link->_faces[0] || !link->_faces[1] )
3253         return debugDumpLink( link );
3254       // a quadrangle sharing <link> with _polygons[iQuad]
3255       _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[iQuad] )];
3256       for ( int i = 0; i < 4; ++i )
3257         if ( quad->_links[i]._link == link )
3258         {
3259           // 1st node of a link opposite to <link> in <quad>
3260           nodes[iL+4] = quad->_links[(i+2)%4].FirstNode();
3261           ++nbN;
3262           break;
3263         }
3264     }
3265     if ( nbN == 8 )
3266       _volumeDefs.set( &nodes[0], 8 );
3267
3268     return nbN == 8;
3269   }
3270   //================================================================================
3271   /*!
3272    * \brief Tries to create a tetrahedron
3273    */
3274   bool Hexahedron::addTetra()
3275   {
3276     int iTria = -1;
3277     for ( size_t i = 0; i < _polygons.size() && iTria < 0; ++i )
3278       if ( _polygons[i]._links.size() == 3 )
3279         iTria = i;
3280     if ( iTria < 0 )
3281       return false;
3282
3283     _Node* nodes[4];
3284     nodes[0] = _polygons[iTria]._links[0].FirstNode();
3285     nodes[1] = _polygons[iTria]._links[1].FirstNode();
3286     nodes[2] = _polygons[iTria]._links[2].FirstNode();
3287
3288     _Link* link = _polygons[iTria]._links[0]._link;
3289     if ( !link->_faces[0] || !link->_faces[1] )
3290       return debugDumpLink( link );
3291
3292     // a triangle sharing <link> with _polygons[0]
3293     _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[iTria] )];
3294     for ( int i = 0; i < 3; ++i )
3295       if ( tria->_links[i]._link == link )
3296       {
3297         nodes[3] = tria->_links[(i+1)%3].LastNode();
3298         _volumeDefs.set( &nodes[0], 4 );
3299         return true;
3300       }
3301
3302     return false;
3303   }
3304   //================================================================================
3305   /*!
3306    * \brief Tries to create a pentahedron
3307    */
3308   bool Hexahedron::addPenta()
3309   {
3310     // find a base triangular face
3311     int iTri = -1;
3312     for ( int iF = 0; iF < 5 && iTri < 0; ++iF )
3313       if ( _polygons[ iF ]._links.size() == 3 )
3314         iTri = iF;
3315     if ( iTri < 0 ) return false;
3316
3317     // find nodes
3318     _Node* nodes[6];
3319     int nbN = 0;
3320     for ( int iL = 0; iL < 3; ++iL )
3321     {
3322       // a base node
3323       nodes[iL] = _polygons[ iTri ]._links[iL].FirstNode();
3324       ++nbN;
3325
3326       // find a top node above the base node
3327       _Link* link = _polygons[ iTri ]._links[iL]._link;
3328       if ( !link->_faces[0] || !link->_faces[1] )
3329         return debugDumpLink( link );
3330       // a quadrangle sharing <link> with a base triangle
3331       _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[ iTri ] )];
3332       if ( quad->_links.size() != 4 ) return false;
3333       for ( int i = 0; i < 4; ++i )
3334         if ( quad->_links[i]._link == link )
3335         {
3336           // 1st node of a link opposite to <link> in <quad>
3337           nodes[iL+3] = quad->_links[(i+2)%4].FirstNode();
3338           ++nbN;
3339           break;
3340         }
3341     }
3342     if ( nbN == 6 )
3343       _volumeDefs.set( &nodes[0], 6 );
3344
3345     return ( nbN == 6 );
3346   }
3347   //================================================================================
3348   /*!
3349    * \brief Tries to create a pyramid
3350    */
3351   bool Hexahedron::addPyra()
3352   {
3353     // find a base quadrangle
3354     int iQuad = -1;
3355     for ( int iF = 0; iF < 5 && iQuad < 0; ++iF )
3356       if ( _polygons[ iF ]._links.size() == 4 )
3357         iQuad = iF;
3358     if ( iQuad < 0 ) return false;
3359
3360     // find nodes
3361     _Node* nodes[5];
3362     nodes[0] = _polygons[iQuad]._links[0].FirstNode();
3363     nodes[1] = _polygons[iQuad]._links[1].FirstNode();
3364     nodes[2] = _polygons[iQuad]._links[2].FirstNode();
3365     nodes[3] = _polygons[iQuad]._links[3].FirstNode();
3366
3367     _Link* link = _polygons[iQuad]._links[0]._link;
3368     if ( !link->_faces[0] || !link->_faces[1] )
3369       return debugDumpLink( link );
3370
3371     // a triangle sharing <link> with a base quadrangle
3372     _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[ iQuad ] )];
3373     if ( tria->_links.size() != 3 ) return false;
3374     for ( int i = 0; i < 3; ++i )
3375       if ( tria->_links[i]._link == link )
3376       {
3377         nodes[4] = tria->_links[(i+1)%3].LastNode();
3378         _volumeDefs.set( &nodes[0], 5 );
3379         return true;
3380       }
3381
3382     return false;
3383   }
3384   //================================================================================
3385   /*!
3386    * \brief Dump a link and return \c false
3387    */
3388   bool Hexahedron::debugDumpLink( Hexahedron::_Link* link )
3389   {
3390 #ifdef _DEBUG_
3391     gp_Pnt p1 = link->_nodes[0]->Point(), p2 = link->_nodes[1]->Point();
3392     cout << "BUG: not shared link. IKJ = ( "<< _i << " " << _j << " " << _k << " )" << endl
3393          << "n1 (" << p1.X() << ", "<< p1.Y() << ", "<< p1.Z() << " )" << endl
3394          << "n2 (" << p2.X() << ", "<< p2.Y() << ", "<< p2.Z() << " )" << endl;
3395 #endif
3396     return false;
3397   }
3398   //================================================================================
3399   /*!
3400    * \brief Classify a point by grid parameters
3401    */
3402   bool Hexahedron::isOutParam(const double uvw[3]) const
3403   {
3404     return (( _grid->_coords[0][ _i   ] - _grid->_tol > uvw[0] ) ||
3405             ( _grid->_coords[0][ _i+1 ] + _grid->_tol < uvw[0] ) ||
3406             ( _grid->_coords[1][ _j   ] - _grid->_tol > uvw[1] ) ||
3407             ( _grid->_coords[1][ _j+1 ] + _grid->_tol < uvw[1] ) ||
3408             ( _grid->_coords[2][ _k   ] - _grid->_tol > uvw[2] ) ||
3409             ( _grid->_coords[2][ _k+1 ] + _grid->_tol < uvw[2] ));
3410   }
3411
3412   //================================================================================
3413   /*!
3414    * \brief computes exact bounding box with axes parallel to given ones
3415    */
3416   //================================================================================
3417
3418   void getExactBndBox( const vector< TopoDS_Shape >& faceVec,
3419                        const double*                 axesDirs,
3420                        Bnd_Box&                      shapeBox )
3421   {
3422     BRep_Builder b;
3423     TopoDS_Compound allFacesComp;
3424     b.MakeCompound( allFacesComp );
3425     for ( size_t iF = 0; iF < faceVec.size(); ++iF )
3426       b.Add( allFacesComp, faceVec[ iF ] );
3427
3428     double sP[6]; // aXmin, aYmin, aZmin, aXmax, aYmax, aZmax
3429     shapeBox.Get(sP[0],sP[1],sP[2],sP[3],sP[4],sP[5]);
3430     double farDist = 0;
3431     for ( int i = 0; i < 6; ++i )
3432       farDist = Max( farDist, 10 * sP[i] );
3433
3434     gp_XYZ axis[3] = { gp_XYZ( axesDirs[0], axesDirs[1], axesDirs[2] ),
3435                        gp_XYZ( axesDirs[3], axesDirs[4], axesDirs[5] ),
3436                        gp_XYZ( axesDirs[6], axesDirs[7], axesDirs[8] ) };
3437     axis[0].Normalize();
3438     axis[1].Normalize();
3439     axis[2].Normalize();
3440
3441     gp_Mat basis( axis[0], axis[1], axis[2] );
3442     gp_Mat bi = basis.Inverted();
3443
3444     gp_Pnt pMin, pMax;
3445     for ( int iDir = 0; iDir < 3; ++iDir )
3446     {
3447       gp_XYZ axis0 = axis[ iDir ];
3448       gp_XYZ axis1 = axis[ ( iDir + 1 ) % 3 ];
3449       gp_XYZ axis2 = axis[ ( iDir + 2 ) % 3 ];
3450       for ( int isMax = 0; isMax < 2; ++isMax )
3451       {
3452         double shift = isMax ? farDist : -farDist;
3453         gp_XYZ orig = shift * axis0;
3454         gp_XYZ norm = axis1 ^ axis2;
3455         gp_Pln pln( orig, norm );
3456         norm = pln.Axis().Direction().XYZ();
3457         BRepBuilderAPI_MakeFace plane( pln, -farDist, farDist, -farDist, farDist );
3458
3459         gp_Pnt& pAxis = isMax ? pMax : pMin;
3460         gp_Pnt pPlane, pFaces;
3461         double dist = GEOMUtils::GetMinDistance( plane, allFacesComp, pPlane, pFaces );
3462         if ( dist < 0 )
3463         {
3464           Bnd_B3d bb;
3465           gp_XYZ corner;
3466           for ( int i = 0; i < 2; ++i ) {
3467             corner.SetCoord( 1, sP[ i*3 ]);
3468             for ( int j = 0; j < 2; ++j ) {
3469               corner.SetCoord( 2, sP[ i*3 + 1 ]);
3470               for ( int k = 0; k < 2; ++k )
3471               {
3472                 corner.SetCoord( 3, sP[ i*3 + 2 ]);
3473                 corner *= bi;
3474                 bb.Add( corner );
3475               }
3476             }
3477           }
3478           corner = isMax ? bb.CornerMax() : bb.CornerMin();
3479           pAxis.SetCoord( iDir+1, corner.Coord( iDir+1 ));
3480         }
3481         else
3482         {
3483           gp_XYZ pf = pFaces.XYZ() * bi;
3484           pAxis.SetCoord( iDir+1, pf.Coord( iDir+1 ) );
3485         }
3486       }
3487     } // loop on 3 axes
3488
3489     shapeBox.SetVoid();
3490     shapeBox.Add( pMin );
3491     shapeBox.Add( pMax );
3492
3493     return;
3494   }
3495
3496 } // namespace
3497
3498 //=============================================================================
3499 /*!
3500  * \brief Generates 3D structured Cartesian mesh in the internal part of
3501  * solid shapes and polyhedral volumes near the shape boundary.
3502  *  \param theMesh - mesh to fill in
3503  *  \param theShape - a compound of all SOLIDs to mesh
3504  *  \retval bool - true in case of success
3505  */
3506 //=============================================================================
3507
3508 bool StdMeshers_Cartesian_3D::Compute(SMESH_Mesh &         theMesh,
3509                                       const TopoDS_Shape & theShape)
3510 {
3511   // The algorithm generates the mesh in following steps:
3512
3513   // 1) Intersection of grid lines with the geometry boundary.
3514   // This step allows to find out if a given node of the initial grid is
3515   // inside or outside the geometry.
3516
3517   // 2) For each cell of the grid, check how many of it's nodes are outside
3518   // of the geometry boundary. Depending on a result of this check
3519   // - skip a cell, if all it's nodes are outside
3520   // - skip a cell, if it is too small according to the size threshold
3521   // - add a hexahedron in the mesh, if all nodes are inside
3522   // - add a polyhedron in the mesh, if some nodes are inside and some outside
3523
3524   _computeCanceled = false;
3525
3526   SMESH_MesherHelper helper( theMesh );
3527
3528   try
3529   {
3530     Grid grid;
3531     grid._helper = &helper;
3532
3533     vector< TopoDS_Shape > faceVec;
3534     {
3535       TopTools_MapOfShape faceMap;
3536       TopExp_Explorer fExp;
3537       for ( fExp.Init( theShape, TopAbs_FACE ); fExp.More(); fExp.Next() )
3538         if ( !faceMap.Add( fExp.Current() ))
3539           faceMap.Remove( fExp.Current() ); // remove a face shared by two solids
3540
3541       for ( fExp.ReInit(); fExp.More(); fExp.Next() )
3542         if ( faceMap.Contains( fExp.Current() ))
3543           faceVec.push_back( fExp.Current() );
3544     }
3545     vector<FaceGridIntersector> facesItersectors( faceVec.size() );
3546     map< TGeomID, vector< TGeomID > > edge2faceIDsMap;
3547     TopExp_Explorer eExp;
3548     Bnd_Box shapeBox;
3549     for ( size_t i = 0; i < faceVec.size(); ++i )
3550     {
3551       facesItersectors[i]._face   = TopoDS::Face    ( faceVec[i] );
3552       facesItersectors[i]._faceID = grid._shapes.Add( faceVec[i] );
3553       facesItersectors[i]._grid   = &grid;
3554       shapeBox.Add( facesItersectors[i].GetFaceBndBox() );
3555
3556       if ( _hyp->GetToAddEdges() )
3557       {
3558         helper.SetSubShape( faceVec[i] );
3559         for ( eExp.Init( faceVec[i], TopAbs_EDGE ); eExp.More(); eExp.Next() )
3560         {
3561           const TopoDS_Edge& edge = TopoDS::Edge( eExp.Current() );
3562           if ( !SMESH_Algo::isDegenerated( edge ) &&
3563                !helper.IsRealSeam( edge ))
3564             edge2faceIDsMap[ grid._shapes.Add( edge )].push_back( facesItersectors[i]._faceID );
3565         }
3566       }
3567     }
3568
3569     getExactBndBox( faceVec, _hyp->GetAxisDirs(), shapeBox );
3570
3571     vector<double> xCoords, yCoords, zCoords;
3572     _hyp->GetCoordinates( xCoords, yCoords, zCoords, shapeBox );
3573
3574     grid.SetCoordinates( xCoords, yCoords, zCoords, _hyp->GetAxisDirs(), shapeBox );
3575
3576     if ( _computeCanceled ) return false;
3577
3578 #ifdef WITH_TBB
3579     { // copy partner faces and curves of not thread-safe types
3580       set< const Standard_Transient* > tshapes;
3581       BRepBuilderAPI_Copy copier;
3582       for ( size_t i = 0; i < facesItersectors.size(); ++i )
3583       {
3584         if ( !facesItersectors[i].IsThreadSafe(tshapes) )
3585         {
3586           copier.Perform( facesItersectors[i]._face );
3587           facesItersectors[i]._face = TopoDS::Face( copier );
3588         }
3589       }
3590     }
3591     // Intersection of grid lines with the geometry boundary.
3592     tbb::parallel_for ( tbb::blocked_range<size_t>( 0, facesItersectors.size() ),
3593                         ParallelIntersector( facesItersectors ),
3594                         tbb::simple_partitioner());
3595 #else
3596     for ( size_t i = 0; i < facesItersectors.size(); ++i )
3597       facesItersectors[i].Intersect();
3598 #endif
3599
3600     // put intersection points onto the GridLine's; this is done after intersection
3601     // to avoid contention of facesItersectors for writing into the same GridLine
3602     // in case of parallel work of facesItersectors
3603     for ( size_t i = 0; i < facesItersectors.size(); ++i )
3604       facesItersectors[i].StoreIntersections();
3605
3606     TopExp_Explorer solidExp (theShape, TopAbs_SOLID);
3607     helper.SetSubShape( solidExp.Current() );
3608     helper.SetElementsOnShape( true );
3609
3610     if ( _computeCanceled ) return false;
3611
3612     // create nodes on the geometry
3613     grid.ComputeNodes(helper);
3614
3615     if ( _computeCanceled ) return false;
3616
3617     // create volume elements
3618     Hexahedron hex( _hyp->GetSizeThreshold(), &grid );
3619     int nbAdded = hex.MakeElements( helper, edge2faceIDsMap );
3620
3621     SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
3622     if ( nbAdded > 0 )
3623     {
3624       // make all SOLIDs computed
3625       if ( SMESHDS_SubMesh* sm1 = meshDS->MeshElements( solidExp.Current()) )
3626       {
3627         SMDS_ElemIteratorPtr volIt = sm1->GetElements();
3628         for ( ; solidExp.More() && volIt->more(); solidExp.Next() )
3629         {
3630           const SMDS_MeshElement* vol = volIt->next();
3631           sm1->RemoveElement( vol );
3632           meshDS->SetMeshElementOnShape( vol, solidExp.Current() );
3633         }
3634       }
3635       // make other sub-shapes computed
3636       setSubmeshesComputed( theMesh, theShape );
3637     }
3638
3639     // remove free nodes
3640     if ( SMESHDS_SubMesh * smDS = meshDS->MeshElements( helper.GetSubShapeID() ))
3641     {
3642       TIDSortedNodeSet nodesToRemove;
3643       // get intersection nodes
3644       for ( int iDir = 0; iDir < 3; ++iDir )
3645       {
3646         vector< GridLine >& lines = grid._lines[ iDir ];
3647         for ( size_t i = 0; i < lines.size(); ++i )
3648         {
3649           multiset< F_IntersectPoint >::iterator ip = lines[i]._intPoints.begin();
3650           for ( ; ip != lines[i]._intPoints.end(); ++ip )
3651             if ( ip->_node && ip->_node->NbInverseElements() == 0 )
3652               nodesToRemove.insert( nodesToRemove.end(), ip->_node );
3653         }
3654       }
3655       // get grid nodes
3656       for ( size_t i = 0; i < grid._nodes.size(); ++i )
3657         if ( grid._nodes[i] && grid._nodes[i]->NbInverseElements() == 0 )
3658           nodesToRemove.insert( nodesToRemove.end(), grid._nodes[i] );
3659
3660       // do remove
3661       TIDSortedNodeSet::iterator n = nodesToRemove.begin();
3662       for ( ; n != nodesToRemove.end(); ++n )
3663         meshDS->RemoveFreeNode( *n, smDS, /*fromGroups=*/false );
3664     }
3665
3666     return nbAdded;
3667
3668   }
3669   // SMESH_ComputeError is not caught at SMESH_submesh level for an unknown reason
3670   catch ( SMESH_ComputeError& e)
3671   {
3672     return error( SMESH_ComputeErrorPtr( new SMESH_ComputeError( e )));
3673   }
3674   return false;
3675 }
3676
3677 //=============================================================================
3678 /*!
3679  *  Evaluate
3680  */
3681 //=============================================================================
3682
3683 bool StdMeshers_Cartesian_3D::Evaluate(SMESH_Mesh &         theMesh,
3684                                        const TopoDS_Shape & theShape,
3685                                        MapShapeNbElems&     theResMap)
3686 {
3687   // TODO
3688 //   std::vector<int> aResVec(SMDSEntity_Last);
3689 //   for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
3690 //   if(IsQuadratic) {
3691 //     aResVec[SMDSEntity_Quad_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
3692 //     int nb1d_face0_int = ( nb2d_face0*4 - nb1d ) / 2;
3693 //     aResVec[SMDSEntity_Node] = nb0d_face0 * ( 2*nb2d/nb1d - 1 ) - nb1d_face0_int * nb2d/nb1d;
3694 //   }
3695 //   else {
3696 //     aResVec[SMDSEntity_Node] = nb0d_face0 * ( nb2d/nb1d - 1 );
3697 //     aResVec[SMDSEntity_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
3698 //   }
3699 //   SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
3700 //   aResMap.insert(std::make_pair(sm,aResVec));
3701
3702   return true;
3703 }
3704
3705 //=============================================================================
3706 namespace
3707 {
3708   /*!
3709    * \brief Event listener setting/unsetting _alwaysComputed flag to
3710    *        submeshes of inferior levels to prevent their computing
3711    */
3712   struct _EventListener : public SMESH_subMeshEventListener
3713   {
3714     string _algoName;
3715
3716     _EventListener(const string& algoName):
3717       SMESH_subMeshEventListener(/*isDeletable=*/true,"StdMeshers_Cartesian_3D::_EventListener"),
3718       _algoName(algoName)
3719     {}
3720     // --------------------------------------------------------------------------------
3721     // setting/unsetting _alwaysComputed flag to submeshes of inferior levels
3722     //
3723     static void setAlwaysComputed( const bool     isComputed,
3724                                    SMESH_subMesh* subMeshOfSolid)
3725     {
3726       SMESH_subMeshIteratorPtr smIt =
3727         subMeshOfSolid->getDependsOnIterator(/*includeSelf=*/false, /*complexShapeFirst=*/false);
3728       while ( smIt->more() )
3729       {
3730         SMESH_subMesh* sm = smIt->next();
3731         sm->SetIsAlwaysComputed( isComputed );
3732       }
3733       subMeshOfSolid->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
3734     }
3735
3736     // --------------------------------------------------------------------------------
3737     // unsetting _alwaysComputed flag if "Cartesian_3D" was removed
3738     //
3739     virtual void ProcessEvent(const int          event,
3740                               const int          eventType,
3741                               SMESH_subMesh*     subMeshOfSolid,
3742                               SMESH_subMeshEventListenerData* data,
3743                               const SMESH_Hypothesis*         hyp = 0)
3744     {
3745       if ( eventType == SMESH_subMesh::COMPUTE_EVENT )
3746       {
3747         setAlwaysComputed( subMeshOfSolid->GetComputeState() == SMESH_subMesh::COMPUTE_OK,
3748                            subMeshOfSolid );
3749       }
3750       else
3751       {
3752         SMESH_Algo* algo3D = subMeshOfSolid->GetAlgo();
3753         if ( !algo3D || _algoName != algo3D->GetName() )
3754           setAlwaysComputed( false, subMeshOfSolid );
3755       }
3756     }
3757
3758     // --------------------------------------------------------------------------------
3759     // set the event listener
3760     //
3761     static void SetOn( SMESH_subMesh* subMeshOfSolid, const string& algoName )
3762     {
3763       subMeshOfSolid->SetEventListener( new _EventListener( algoName ),
3764                                         /*data=*/0,
3765                                         subMeshOfSolid );
3766     }
3767
3768   }; // struct _EventListener
3769
3770 } // namespace
3771
3772 //================================================================================
3773 /*!
3774  * \brief Sets event listener to submeshes if necessary
3775  *  \param subMesh - submesh where algo is set
3776  * This method is called when a submesh gets HYP_OK algo_state.
3777  * After being set, event listener is notified on each event of a submesh.
3778  */
3779 //================================================================================
3780
3781 void StdMeshers_Cartesian_3D::SetEventListener(SMESH_subMesh* subMesh)
3782 {
3783   _EventListener::SetOn( subMesh, GetName() );
3784 }
3785
3786 //================================================================================
3787 /*!
3788  * \brief Set _alwaysComputed flag to submeshes of inferior levels to avoid their computing
3789  */
3790 //================================================================================
3791
3792 void StdMeshers_Cartesian_3D::setSubmeshesComputed(SMESH_Mesh&         theMesh,
3793                                                    const TopoDS_Shape& theShape)
3794 {
3795   for ( TopExp_Explorer soExp( theShape, TopAbs_SOLID ); soExp.More(); soExp.Next() )
3796     _EventListener::setAlwaysComputed( true, theMesh.GetSubMesh( soExp.Current() ));
3797 }
3798