Salome HOME
#17237: Body fitting on sub-mesh, #16523: Treatment of internal faces
[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 #include "StdMeshers_CartesianParameters3D.hxx"
27
28 #include "ObjectPool.hxx"
29 #include "SMDS_MeshNode.hxx"
30 #include "SMDS_VolumeTool.hxx"
31 #include "SMESHDS_Mesh.hxx"
32 #include "SMESH_Block.hxx"
33 #include "SMESH_Comment.hxx"
34 #include "SMESH_ControlsDef.hxx"
35 #include "SMESH_Mesh.hxx"
36 #include "SMESH_MeshAlgos.hxx"
37 #include "SMESH_MeshEditor.hxx"
38 #include "SMESH_MesherHelper.hxx"
39 #include "SMESH_subMesh.hxx"
40 #include "SMESH_subMeshEventListener.hxx"
41 #include "StdMeshers_FaceSide.hxx"
42
43 #include <utilities.h>
44 #include <Utils_ExceptHandlers.hxx>
45
46 #include <GEOMUtils.hxx>
47
48 #include <BRepAdaptor_Curve.hxx>
49 #include <BRepAdaptor_Surface.hxx>
50 #include <BRepBndLib.hxx>
51 #include <BRepBuilderAPI_Copy.hxx>
52 #include <BRepBuilderAPI_MakeFace.hxx>
53 #include <BRepTools.hxx>
54 #include <BRepTopAdaptor_FClass2d.hxx>
55 #include <BRep_Builder.hxx>
56 #include <BRep_Tool.hxx>
57 #include <Bnd_B3d.hxx>
58 #include <Bnd_Box.hxx>
59 #include <ElSLib.hxx>
60 #include <GCPnts_UniformDeflection.hxx>
61 #include <Geom2d_BSplineCurve.hxx>
62 #include <Geom2d_BezierCurve.hxx>
63 #include <Geom2d_TrimmedCurve.hxx>
64 #include <GeomAPI_ProjectPointOnSurf.hxx>
65 #include <GeomLib.hxx>
66 #include <Geom_BSplineCurve.hxx>
67 #include <Geom_BSplineSurface.hxx>
68 #include <Geom_BezierCurve.hxx>
69 #include <Geom_BezierSurface.hxx>
70 #include <Geom_RectangularTrimmedSurface.hxx>
71 #include <Geom_TrimmedCurve.hxx>
72 #include <IntAna_IntConicQuad.hxx>
73 #include <IntAna_IntLinTorus.hxx>
74 #include <IntAna_Quadric.hxx>
75 #include <IntCurveSurface_TransitionOnCurve.hxx>
76 #include <IntCurvesFace_Intersector.hxx>
77 #include <Poly_Triangulation.hxx>
78 #include <Precision.hxx>
79 #include <TopExp.hxx>
80 #include <TopExp_Explorer.hxx>
81 #include <TopLoc_Location.hxx>
82 #include <TopTools_IndexedMapOfShape.hxx>
83 #include <TopTools_MapOfShape.hxx>
84 #include <TopoDS.hxx>
85 #include <TopoDS_Compound.hxx>
86 #include <TopoDS_Face.hxx>
87 #include <TopoDS_TShape.hxx>
88 #include <gp_Cone.hxx>
89 #include <gp_Cylinder.hxx>
90 #include <gp_Lin.hxx>
91 #include <gp_Pln.hxx>
92 #include <gp_Pnt2d.hxx>
93 #include <gp_Sphere.hxx>
94 #include <gp_Torus.hxx>
95
96 #include <limits>
97
98 #include <boost/container/flat_map.hpp>
99
100 //#undef WITH_TBB
101 #ifdef WITH_TBB
102
103 #ifdef WIN32
104 // See https://docs.microsoft.com/en-gb/cpp/porting/modifying-winver-and-win32-winnt?view=vs-2019
105 // Windows 10 = 0x0A00  
106 #define WINVER 0x0A00
107 #define _WIN32_WINNT 0x0A00
108 #endif
109
110 #include <tbb/parallel_for.h>
111 //#include <tbb/enumerable_thread_specific.h>
112 #endif
113
114 using namespace std;
115 using namespace SMESH;
116
117 #ifdef _DEBUG_
118 //#define _MY_DEBUG_
119 #endif
120
121 //=============================================================================
122 /*!
123  * Constructor
124  */
125 //=============================================================================
126
127 StdMeshers_Cartesian_3D::StdMeshers_Cartesian_3D(int hypId, SMESH_Gen * gen)
128   :SMESH_3D_Algo(hypId, gen)
129 {
130   _name = "Cartesian_3D";
131   _shapeType = (1 << TopAbs_SOLID);       // 1 bit /shape type
132   _compatibleHypothesis.push_back("CartesianParameters3D");
133
134   _onlyUnaryInput = false;          // to mesh all SOLIDs at once
135   _requireDiscreteBoundary = false; // 2D mesh not needed
136   _supportSubmeshes = false;        // do not use any existing mesh
137 }
138
139 //=============================================================================
140 /*!
141  * Check presence of a hypothesis
142  */
143 //=============================================================================
144
145 bool StdMeshers_Cartesian_3D::CheckHypothesis (SMESH_Mesh&          aMesh,
146                                                const TopoDS_Shape&  aShape,
147                                                Hypothesis_Status&   aStatus)
148 {
149   aStatus = SMESH_Hypothesis::HYP_MISSING;
150
151   const list<const SMESHDS_Hypothesis*>& hyps = GetUsedHypothesis(aMesh, aShape);
152   list <const SMESHDS_Hypothesis* >::const_iterator h = hyps.begin();
153   if ( h == hyps.end())
154   {
155     return false;
156   }
157
158   for ( ; h != hyps.end(); ++h )
159   {
160     if (( _hyp = dynamic_cast<const StdMeshers_CartesianParameters3D*>( *h )))
161     {
162       aStatus = _hyp->IsDefined() ? HYP_OK : HYP_BAD_PARAMETER;
163       break;
164     }
165   }
166
167   return aStatus == HYP_OK;
168 }
169
170 namespace
171 {
172   typedef int TGeomID; // IDs of sub-shapes
173
174   //=============================================================================
175   // Definitions of internal utils
176   // --------------------------------------------------------------------------
177   enum Transition {
178     Trans_TANGENT = IntCurveSurface_Tangent,
179     Trans_IN      = IntCurveSurface_In,
180     Trans_OUT     = IntCurveSurface_Out,
181     Trans_APEX,
182     Trans_INTERNAL // for INTERNAL FACE
183   };
184   // --------------------------------------------------------------------------
185   /*!
186    * \brief Container of IDs of SOLID sub-shapes
187    */
188   class Solid // sole SOLID contains all sub-shapes
189   {
190     TGeomID _id; // SOLID id
191     bool    _hasInternalFaces;
192   public:
193     virtual ~Solid() {}
194     virtual bool Contains( TGeomID subID ) const { return true; }
195     virtual bool ContainsAny( const vector< TGeomID>& subIDs ) const { return true; }
196     virtual TopAbs_Orientation Orientation( const TopoDS_Shape& s ) const { return s.Orientation(); }
197     virtual bool IsOutsideOriented( TGeomID faceID ) const { return true; }
198     void SetID( TGeomID id ) { _id = id; }
199     TGeomID ID() const { return _id; }
200     void SetHasInternalFaces( bool has ) { _hasInternalFaces = has; }
201     bool HasInternalFaces() const { return _hasInternalFaces; }
202   };
203   // --------------------------------------------------------------------------
204   class OneOfSolids : public Solid
205   {
206     TColStd_MapOfInteger _subIDs;
207     TopTools_MapOfShape  _faces; // keep FACE orientation
208     TColStd_MapOfInteger _outFaceIDs; // FACEs of shape_to_mesh oriented outside the SOLID
209   public:
210     void Init( const TopoDS_Shape& solid,
211                TopAbs_ShapeEnum    subType,
212                const SMESHDS_Mesh* mesh );
213     virtual bool Contains( TGeomID i ) const { return i == ID() || _subIDs.Contains( i ); }
214     virtual bool ContainsAny( const vector< TGeomID>& subIDs ) const
215     {
216       for ( size_t i = 0; i < subIDs.size(); ++i ) if ( Contains( subIDs[ i ])) return true;
217       return false;
218     }
219     virtual TopAbs_Orientation Orientation( const TopoDS_Shape& face ) const
220     {
221       const TopoDS_Shape& sInMap = const_cast< OneOfSolids* >(this)->_faces.Added( face );
222       return sInMap.Orientation();
223     }
224     virtual bool IsOutsideOriented( TGeomID faceID ) const
225     {
226       return faceID == 0 || _outFaceIDs.Contains( faceID );
227     }
228   };
229   // --------------------------------------------------------------------------
230   /*!
231    * \brief Geom data
232    */
233   struct Geometry
234   {
235     TopoDS_Shape                _mainShape;
236     vector< vector< TGeomID > > _solidIDsByShapeID;// V/E/F ID -> SOLID IDs
237     Solid                       _soleSolid;
238     map< TGeomID, OneOfSolids > _solidByID;
239     TColStd_MapOfInteger        _boundaryFaces; // FACEs on boundary of mesh->ShapeToMesh()
240     TColStd_MapOfInteger        _strangeEdges; // EDGEs shared by strange FACEs
241     TGeomID                     _extIntFaceID; // pseudo FACE - extension of INTERNAL FACE
242
243     Controls::ElementsOnShape _edgeClassifier;
244     Controls::ElementsOnShape _vertexClassifier;
245
246     bool IsOneSolid() const { return _solidByID.size() < 2; }
247   };
248   // --------------------------------------------------------------------------
249   /*!
250    * \brief Common data of any intersection between a Grid and a shape
251    */
252   struct B_IntersectPoint
253   {
254     mutable const SMDS_MeshNode* _node;
255     mutable vector< TGeomID >    _faceIDs;
256
257     B_IntersectPoint(): _node(NULL) {}
258     void Add( const vector< TGeomID >& fIDs, const SMDS_MeshNode* n=0 ) const;
259     int HasCommonFace( const B_IntersectPoint * other, int avoidFace=-1 ) const;
260     bool IsOnFace( int faceID ) const;
261     virtual ~B_IntersectPoint() {}
262   };
263   // --------------------------------------------------------------------------
264   /*!
265    * \brief Data of intersection between a GridLine and a TopoDS_Face
266    */
267   struct F_IntersectPoint : public B_IntersectPoint
268   {
269     double             _paramOnLine;
270     double             _u, _v;
271     mutable Transition _transition;
272     mutable size_t     _indexOnLine;
273
274     bool operator< ( const F_IntersectPoint& o ) const { return _paramOnLine < o._paramOnLine; }
275   };
276   // --------------------------------------------------------------------------
277   /*!
278    * \brief Data of intersection between GridPlanes and a TopoDS_EDGE
279    */
280   struct E_IntersectPoint : public B_IntersectPoint
281   {
282     gp_Pnt  _point;
283     double  _uvw[3];
284     TGeomID _shapeID; // ID of EDGE or VERTEX
285   };
286   // --------------------------------------------------------------------------
287   /*!
288    * \brief A line of the grid and its intersections with 2D geometry
289    */
290   struct GridLine
291   {
292     gp_Lin _line;
293     double _length; // line length
294     multiset< F_IntersectPoint > _intPoints;
295
296     void RemoveExcessIntPoints( const double tol );
297     TGeomID GetSolidIDBefore( multiset< F_IntersectPoint >::iterator ip,
298                               const TGeomID                          prevID,
299                               const Geometry&                        geom);
300   };
301   // --------------------------------------------------------------------------
302   /*!
303    * \brief Planes of the grid used to find intersections of an EDGE with a hexahedron
304    */
305   struct GridPlanes
306   {
307     gp_XYZ           _zNorm;
308     vector< gp_XYZ > _origins; // origin points of all planes in one direction
309     vector< double > _zProjs;  // projections of origins to _zNorm
310   };
311   // --------------------------------------------------------------------------
312   /*!
313    * \brief Iterator on the parallel grid lines of one direction
314    */
315   struct LineIndexer
316   {
317     size_t _size  [3];
318     size_t _curInd[3];
319     size_t _iVar1, _iVar2, _iConst;
320     string _name1, _name2, _nameConst;
321     LineIndexer() {}
322     LineIndexer( size_t sz1, size_t sz2, size_t sz3,
323                  size_t iv1, size_t iv2, size_t iConst,
324                  const string& nv1, const string& nv2, const string& nConst )
325     {
326       _size[0] = sz1; _size[1] = sz2; _size[2] = sz3;
327       _curInd[0] = _curInd[1] = _curInd[2] = 0;
328       _iVar1 = iv1; _iVar2 = iv2; _iConst = iConst;
329       _name1 = nv1; _name2 = nv2; _nameConst = nConst;
330     }
331
332     size_t I() const { return _curInd[0]; }
333     size_t J() const { return _curInd[1]; }
334     size_t K() const { return _curInd[2]; }
335     void SetIJK( size_t i, size_t j, size_t k )
336     {
337       _curInd[0] = i; _curInd[1] = j; _curInd[2] = k;
338     }
339     void operator++()
340     {
341       if ( ++_curInd[_iVar1] == _size[_iVar1] )
342         _curInd[_iVar1] = 0, ++_curInd[_iVar2];
343     }
344     bool More() const { return _curInd[_iVar2] < _size[_iVar2]; }
345     size_t LineIndex   () const { return _curInd[_iVar1] + _curInd[_iVar2]* _size[_iVar1]; }
346     size_t LineIndex10 () const { return (_curInd[_iVar1] + 1 ) + _curInd[_iVar2]* _size[_iVar1]; }
347     size_t LineIndex01 () const { return _curInd[_iVar1] + (_curInd[_iVar2] + 1 )* _size[_iVar1]; }
348     size_t LineIndex11 () const { return (_curInd[_iVar1] + 1 ) + (_curInd[_iVar2] + 1 )* _size[_iVar1]; }
349     void SetIndexOnLine (size_t i)  { _curInd[ _iConst ] = i; }
350     size_t NbLines() const { return _size[_iVar1] * _size[_iVar2]; }
351   };
352   // --------------------------------------------------------------------------
353   /*!
354    * \brief Container of GridLine's
355    */
356   struct Grid
357   {
358     vector< double >   _coords[3]; // coordinates of grid nodes
359     gp_XYZ             _axes  [3]; // axis directions
360     vector< GridLine > _lines [3]; //    in 3 directions
361     double             _tol, _minCellSize;
362     gp_XYZ             _origin;
363     gp_Mat             _invB; // inverted basis of _axes
364
365     vector< const SMDS_MeshNode* >    _nodes; // mesh nodes at grid nodes
366     vector< const F_IntersectPoint* > _gridIntP; // grid node intersection with geometry
367     ObjectPool< E_IntersectPoint >    _edgeIntPool; // intersections with EDGEs
368     ObjectPool< F_IntersectPoint >    _extIntPool; // intersections with extended INTERNAL FACEs
369     //list< E_IntersectPoint >          _edgeIntP; // intersections with EDGEs
370
371     Geometry                          _geometry;
372     bool                              _toAddEdges;
373     bool                              _toCreateFaces;
374     bool                              _toConsiderInternalFaces;
375     bool                              _toUseThresholdForInternalFaces;
376     double                            _sizeThreshold;
377
378     SMESH_MesherHelper*               _helper;
379
380     size_t CellIndex( size_t i, size_t j, size_t k ) const
381     {
382       return i + j*(_coords[0].size()-1) + k*(_coords[0].size()-1)*(_coords[1].size()-1);
383     }
384     size_t NodeIndex( size_t i, size_t j, size_t k ) const
385     {
386       return i + j*_coords[0].size() + k*_coords[0].size()*_coords[1].size();
387     }
388     size_t NodeIndexDX() const { return 1; }
389     size_t NodeIndexDY() const { return _coords[0].size(); }
390     size_t NodeIndexDZ() const { return _coords[0].size() * _coords[1].size(); }
391
392     LineIndexer GetLineIndexer(size_t iDir) const;
393
394     E_IntersectPoint* Add( const E_IntersectPoint& ip )
395     {
396       E_IntersectPoint* eip = _edgeIntPool.getNew();
397       *eip = ip;
398       return eip;
399     }
400     void Remove( E_IntersectPoint* eip ) { _edgeIntPool.destroy( eip ); }
401
402     TGeomID ShapeID( const TopoDS_Shape& s ) const;
403     const TopoDS_Shape& Shape( TGeomID id ) const;
404     TopAbs_ShapeEnum ShapeType( TGeomID id ) const { return Shape(id).ShapeType(); }
405     void InitGeometry( const TopoDS_Shape& theShape );
406     void InitClassifier( const TopoDS_Shape&        mainShape,
407                          TopAbs_ShapeEnum           shapeType,
408                          Controls::ElementsOnShape& classifier );
409     void GetEdgesToImplement( map< TGeomID, vector< TGeomID > > & edge2faceMap,
410                               const TopoDS_Shape&                 shape,
411                               const vector< TopoDS_Shape >&       faces );
412     void SetSolidFather( const TopoDS_Shape& s, const TopoDS_Shape& theShapeToMesh );
413     bool IsShared( TGeomID faceID ) const;
414     bool IsAnyShared( const std::vector< TGeomID >& faceIDs ) const;
415     bool IsInternal( TGeomID faceID ) const {
416       return ( faceID == PseudoIntExtFaceID() ||
417                Shape( faceID ).Orientation() == TopAbs_INTERNAL ); }
418     bool IsSolid( TGeomID shapeID ) const {
419       if ( _geometry.IsOneSolid() ) return _geometry._soleSolid.ID() == shapeID;
420       else                          return _geometry._solidByID.count( shapeID ); }
421     bool IsStrangeEdge( TGeomID id ) const { return _geometry._strangeEdges.Contains( id ); }
422     TGeomID PseudoIntExtFaceID() const { return _geometry._extIntFaceID; }
423     Solid* GetSolid( TGeomID solidID = 0 );
424     Solid* GetOneOfSolids( TGeomID solidID );
425     const vector< TGeomID > & GetSolidIDs( TGeomID subShapeID ) const;
426     bool IsCorrectTransition( TGeomID faceID, const Solid* solid );
427     bool IsBoundaryFace( TGeomID face ) const { return _geometry._boundaryFaces.Contains( face ); }
428     void SetOnShape( const SMDS_MeshNode* n, const F_IntersectPoint& ip, bool unset=false );
429     bool IsToCheckNodePos() const { return !_toAddEdges && _toCreateFaces; }
430
431     void SetCoordinates(const vector<double>& xCoords,
432                         const vector<double>& yCoords,
433                         const vector<double>& zCoords,
434                         const double*         axesDirs,
435                         const Bnd_Box&        bndBox );
436     void ComputeUVW(const gp_XYZ& p, double uvw[3]);
437     void ComputeNodes(SMESH_MesherHelper& helper);
438   };
439   // --------------------------------------------------------------------------
440   /*!
441    * \brief Return cells sharing a link
442    */
443   struct CellsAroundLink
444   {
445     int    _dInd[4][3];
446     size_t _nbCells[3];
447     int    _i,_j,_k;
448     Grid*  _grid;
449
450     CellsAroundLink( Grid* grid, int iDir ):
451       _dInd{ {0,0,0}, {0,0,0}, {0,0,0}, {0,0,0} },
452       _nbCells{ grid->_coords[0].size() - 1,
453           grid->_coords[1].size() - 1,
454           grid->_coords[2].size() - 1 },
455       _grid( grid )
456     {
457       const int iDirOther[3][2] = {{ 1,2 },{ 0,2 },{ 0,1 }};
458       _dInd[1][ iDirOther[iDir][0] ] = -1;
459       _dInd[2][ iDirOther[iDir][1] ] = -1;
460       _dInd[3][ iDirOther[iDir][0] ] = -1; _dInd[3][ iDirOther[iDir][1] ] = -1;
461     }
462     void Init( int i, int j, int k, int link12 = 0 )
463     {
464       int iL = link12 % 4;
465       _i = i - _dInd[iL][0];
466       _j = j - _dInd[iL][1];
467       _k = k - _dInd[iL][2];
468     }
469     bool GetCell( int iL, int& i, int& j, int& k, int& cellIndex )
470     {
471       i =  _i + _dInd[iL][0];
472       j =  _j + _dInd[iL][1];
473       k =  _k + _dInd[iL][2];
474       if ( i < 0 || i >= (int)_nbCells[0] ||
475            j < 0 || j >= (int)_nbCells[1] ||
476            k < 0 || k >= (int)_nbCells[2] )
477         return false;
478       cellIndex = _grid->CellIndex( i,j,k );
479       return true;
480     }
481   };
482   // --------------------------------------------------------------------------
483   /*!
484    * \brief Intersector of TopoDS_Face with all GridLine's
485    */
486   struct FaceGridIntersector
487   {
488     TopoDS_Face _face;
489     TGeomID     _faceID;
490     Grid*       _grid;
491     Bnd_Box     _bndBox;
492     IntCurvesFace_Intersector* _surfaceInt;
493     vector< std::pair< GridLine*, F_IntersectPoint > > _intersections;
494
495     FaceGridIntersector(): _grid(0), _surfaceInt(0) {}
496     void Intersect();
497
498     void StoreIntersections()
499     {
500       for ( size_t i = 0; i < _intersections.size(); ++i )
501       {
502         multiset< F_IntersectPoint >::iterator ip =
503           _intersections[i].first->_intPoints.insert( _intersections[i].second );
504         ip->_faceIDs.reserve( 1 );
505         ip->_faceIDs.push_back( _faceID );
506       }
507     }
508     const Bnd_Box& GetFaceBndBox()
509     {
510       GetCurveFaceIntersector();
511       return _bndBox;
512     }
513     IntCurvesFace_Intersector* GetCurveFaceIntersector()
514     {
515       if ( !_surfaceInt )
516       {
517         _surfaceInt = new IntCurvesFace_Intersector( _face, Precision::PConfusion() );
518         _bndBox     = _surfaceInt->Bounding();
519         if ( _bndBox.IsVoid() )
520           BRepBndLib::Add (_face, _bndBox);
521       }
522       return _surfaceInt;
523     }
524     bool IsThreadSafe(set< const Standard_Transient* >& noSafeTShapes) const;
525   };
526   // --------------------------------------------------------------------------
527   /*!
528    * \brief Intersector of a surface with a GridLine
529    */
530   struct FaceLineIntersector
531   {
532     double      _tol;
533     double      _u, _v, _w; // params on the face and the line
534     Transition  _transition; // transition at intersection (see IntCurveSurface.cdl)
535     Transition  _transIn, _transOut; // IN and OUT transitions depending of face orientation
536
537     gp_Pln      _plane;
538     gp_Cylinder _cylinder;
539     gp_Cone     _cone;
540     gp_Sphere   _sphere;
541     gp_Torus    _torus;
542     IntCurvesFace_Intersector* _surfaceInt;
543
544     vector< F_IntersectPoint > _intPoints;
545
546     void IntersectWithPlane   (const GridLine& gridLine);
547     void IntersectWithCylinder(const GridLine& gridLine);
548     void IntersectWithCone    (const GridLine& gridLine);
549     void IntersectWithSphere  (const GridLine& gridLine);
550     void IntersectWithTorus   (const GridLine& gridLine);
551     void IntersectWithSurface (const GridLine& gridLine);
552
553     bool UVIsOnFace() const;
554     void addIntPoint(const bool toClassify=true);
555     bool isParamOnLineOK( const double linLength )
556     {
557       return -_tol < _w && _w < linLength + _tol;
558     }
559     FaceLineIntersector():_surfaceInt(0) {}
560     ~FaceLineIntersector() { if (_surfaceInt ) delete _surfaceInt; _surfaceInt = 0; }
561   };
562   // --------------------------------------------------------------------------
563   /*!
564    * \brief Class representing topology of the hexahedron and creating a mesh
565    *        volume basing on analysis of hexahedron intersection with geometry
566    */
567   class Hexahedron
568   {
569     // --------------------------------------------------------------------------------
570     struct _Face;
571     struct _Link;
572     enum IsInternalFlag { IS_NOT_INTERNAL, IS_INTERNAL, IS_CUT_BY_INTERNAL_FACE };
573     // --------------------------------------------------------------------------------
574     struct _Node //!< node either at a hexahedron corner or at intersection
575     {
576       const SMDS_MeshNode*    _node; // mesh node at hexahedron corner
577       const B_IntersectPoint* _intPoint;
578       const _Face*            _usedInFace;
579       char                    _isInternalFlags;
580
581       _Node(const SMDS_MeshNode* n=0, const B_IntersectPoint* ip=0)
582         :_node(n), _intPoint(ip), _usedInFace(0), _isInternalFlags(0) {} 
583       const SMDS_MeshNode*    Node() const
584       { return ( _intPoint && _intPoint->_node ) ? _intPoint->_node : _node; }
585       const E_IntersectPoint* EdgeIntPnt() const
586       { return static_cast< const E_IntersectPoint* >( _intPoint ); }
587       const F_IntersectPoint* FaceIntPnt() const
588       { return static_cast< const F_IntersectPoint* >( _intPoint ); }
589       const vector< TGeomID >& faces() const { return _intPoint->_faceIDs; }
590       TGeomID face(size_t i) const { return _intPoint->_faceIDs[ i ]; }
591       void SetInternal( IsInternalFlag intFlag ) { _isInternalFlags |= intFlag; }
592       bool IsCutByInternal() const { return _isInternalFlags & IS_CUT_BY_INTERNAL_FACE; }
593       bool IsUsedInFace( const _Face* polygon = 0 )
594       {
595         return polygon ? ( _usedInFace == polygon ) : bool( _usedInFace );
596       }
597       TGeomID IsLinked( const B_IntersectPoint* other,
598                         TGeomID                 avoidFace=-1 ) const // returns id of a common face
599       {
600         return _intPoint ? _intPoint->HasCommonFace( other, avoidFace ) : 0;
601       }
602       bool IsOnFace( TGeomID faceID ) const // returns true if faceID is found
603       {
604         return _intPoint ? _intPoint->IsOnFace( faceID ) : false;
605       }
606       gp_Pnt Point() const
607       {
608         if ( const SMDS_MeshNode* n = Node() )
609           return SMESH_NodeXYZ( n );
610         if ( const E_IntersectPoint* eip =
611              dynamic_cast< const E_IntersectPoint* >( _intPoint ))
612           return eip->_point;
613         return gp_Pnt( 1e100, 0, 0 );
614       }
615       TGeomID ShapeID() const
616       {
617         if ( const E_IntersectPoint* eip = dynamic_cast< const E_IntersectPoint* >( _intPoint ))
618           return eip->_shapeID;
619         return 0;
620       }
621       void Add( const E_IntersectPoint* ip )
622       {
623         // Possible cases before Add(ip):
624         ///  1) _node != 0 --> _Node at hex corner ( _intPoint == 0 || _intPoint._node == 0 )
625         ///  2) _node == 0 && _intPoint._node != 0  -->  link intersected by FACE
626         ///  3) _node == 0 && _intPoint._node == 0  -->  _Node at EDGE intersection
627         //
628         // If ip is added in cases 1) and 2) _node position must be changed to ip._shapeID
629         //   at creation of elements
630         // To recognize this case, set _intPoint._node = Node()
631         const SMDS_MeshNode* node = Node();
632         if ( !_intPoint ) {
633           _intPoint = ip;
634         }
635         else {
636           ip->Add( _intPoint->_faceIDs );
637           _intPoint = ip;
638         }
639         if ( node )
640           _node = _intPoint->_node = node;
641       }
642     };
643     // --------------------------------------------------------------------------------
644     struct _Link // link connecting two _Node's
645     {
646       _Node* _nodes[2];
647       _Face* _faces[2]; // polygons sharing a link
648       vector< const F_IntersectPoint* > _fIntPoints; // GridLine intersections with FACEs
649       vector< _Node* >                  _fIntNodes;   // _Node's at _fIntPoints
650       vector< _Link >                   _splits;
651       _Link(): _faces{ 0, 0 } {}
652     };
653     // --------------------------------------------------------------------------------
654     struct _OrientedLink
655     {
656       _Link* _link;
657       bool   _reverse;
658       _OrientedLink( _Link* link=0, bool reverse=false ): _link(link), _reverse(reverse) {}
659       void Reverse() { _reverse = !_reverse; }
660       int NbResultLinks() const { return _link->_splits.size(); }
661       _OrientedLink ResultLink(int i) const
662       {
663         return _OrientedLink(&_link->_splits[_reverse ? NbResultLinks()-i-1 : i],_reverse);
664       }
665       _Node* FirstNode() const { return _link->_nodes[ _reverse ]; }
666       _Node* LastNode()  const { return _link->_nodes[ !_reverse ]; }
667       operator bool() const { return _link; }
668       vector< TGeomID > GetNotUsedFace(const set<TGeomID>& usedIDs ) const // returns supporting FACEs
669       {
670         vector< TGeomID > faces;
671         const B_IntersectPoint *ip0, *ip1;
672         if (( ip0 = _link->_nodes[0]->_intPoint ) &&
673             ( ip1 = _link->_nodes[1]->_intPoint ))
674         {
675           for ( size_t i = 0; i < ip0->_faceIDs.size(); ++i )
676             if ( ip1->IsOnFace ( ip0->_faceIDs[i] ) &&
677                  !usedIDs.count( ip0->_faceIDs[i] ) )
678               faces.push_back( ip0->_faceIDs[i] );
679         }
680         return faces;
681       }
682       bool HasEdgeNodes() const
683       {
684         return ( dynamic_cast< const E_IntersectPoint* >( _link->_nodes[0]->_intPoint ) ||
685                  dynamic_cast< const E_IntersectPoint* >( _link->_nodes[1]->_intPoint ));
686       }
687       int NbFaces() const
688       {
689         return !_link->_faces[0] ? 0 : 1 + bool( _link->_faces[1] );
690       }
691       void AddFace( _Face* f )
692       {
693         if ( _link->_faces[0] )
694         {
695           _link->_faces[1] = f;
696         }
697         else
698         {
699           _link->_faces[0] = f;
700           _link->_faces[1] = 0;
701         }
702       }
703       void RemoveFace( _Face* f )
704       {
705         if ( !_link->_faces[0] ) return;
706
707         if ( _link->_faces[1] == f )
708         {
709           _link->_faces[1] = 0;
710         }
711         else if ( _link->_faces[0] == f )
712         {
713           _link->_faces[0] = 0;
714           if ( _link->_faces[1] )
715           {
716             _link->_faces[0] = _link->_faces[1];
717             _link->_faces[1] = 0;
718           }
719         }
720       }
721     };
722     // --------------------------------------------------------------------------------
723     struct _SplitIterator //! set to _hexLinks splits on one side of INTERNAL FACEs
724     {
725       struct _Split // data of a link split
726       {
727         int    _linkID;      // hex link ID
728         _Node* _nodes[2];
729         int    _iCheckIteration; // iteration where split is tried as Hexahedron split
730         _Link* _checkedSplit;    // split set to hex links
731         bool   _isUsed;    // used in a volume
732
733         _Split( _Link & split, int iLink ):
734           _linkID( iLink ), _nodes{ split._nodes[0], split._nodes[1] },
735           _iCheckIteration( 0 ), _isUsed( false )
736         {}
737         bool IsCheckedOrUsed( bool used ) const { return used ? _isUsed : _iCheckIteration > 0; }
738       };
739       _Link*                _hexLinks;
740       std::vector< _Split > _splits;
741       int                   _iterationNb;
742       size_t                _nbChecked;
743       size_t                _nbUsed;
744       std::vector< _Node* > _freeNodes; // nodes reached while composing a split set
745
746       _SplitIterator( _Link* hexLinks ):
747         _hexLinks( hexLinks ), _iterationNb(0), _nbChecked(0), _nbUsed(0)
748       {
749         _freeNodes.reserve( 12 );
750         _splits.reserve( 24 );
751         for ( int iL = 0; iL < 12; ++iL )
752           for ( size_t iS = 0; iS < _hexLinks[ iL ]._splits.size(); ++iS )
753             _splits.emplace_back( _hexLinks[ iL ]._splits[ iS ], iL );
754         Next();
755       }
756       bool More() const { return _nbUsed < _splits.size(); }
757       bool Next();
758     };
759     // --------------------------------------------------------------------------------
760     struct _Face
761     {
762       vector< _OrientedLink > _links;       // links on GridLine's
763       vector< _Link >         _polyLinks;   // links added to close a polygonal face
764       vector< _Node* >        _eIntNodes;   // nodes at intersection with EDGEs
765       bool IsPolyLink( const _OrientedLink& ol )
766       {
767         return _polyLinks.empty() ? false :
768           ( &_polyLinks[0] <= ol._link &&  ol._link <= &_polyLinks.back() );
769       }
770       void AddPolyLink(_Node* n0, _Node* n1, _Face* faceToFindEqual=0)
771       {
772         if ( faceToFindEqual && faceToFindEqual != this ) {
773           for ( size_t iL = 0; iL < faceToFindEqual->_polyLinks.size(); ++iL )
774             if ( faceToFindEqual->_polyLinks[iL]._nodes[0] == n1 &&
775                  faceToFindEqual->_polyLinks[iL]._nodes[1] == n0 )
776             {
777               _links.push_back
778                 ( _OrientedLink( & faceToFindEqual->_polyLinks[iL], /*reverse=*/true ));
779               return;
780             }
781         }
782         _Link l;
783         l._nodes[0] = n0;
784         l._nodes[1] = n1;
785         _polyLinks.push_back( l );
786         _links.push_back( _OrientedLink( &_polyLinks.back() ));
787       }
788     };
789     // --------------------------------------------------------------------------------
790     struct _volumeDef // holder of nodes of a volume mesh element
791     {
792       struct _nodeDef
793       {
794         const SMDS_MeshNode*    _node; // mesh node at hexahedron corner
795         const B_IntersectPoint* _intPoint;
796
797         _nodeDef( _Node* n ): _node( n->_node), _intPoint( n->_intPoint ) {}
798         const SMDS_MeshNode*    Node() const
799         { return ( _intPoint && _intPoint->_node ) ? _intPoint->_node : _node; }
800         const E_IntersectPoint* EdgeIntPnt() const
801         { return static_cast< const E_IntersectPoint* >( _intPoint ); }
802       };
803       vector< _nodeDef >      _nodes;
804       vector< int >           _quantities;
805       _volumeDef*             _next; // to store several _volumeDefs in a chain
806       TGeomID                 _solidID;
807       const SMDS_MeshElement* _volume; // new volume
808
809       _volumeDef(): _next(0), _solidID(0), _volume(0) {}
810       ~_volumeDef() { delete _next; }
811       _volumeDef( _volumeDef& other ):
812         _next(0), _solidID( other._solidID ), _volume( other._volume )
813       { _nodes.swap( other._nodes ); _quantities.swap( other._quantities ); other._volume = 0; }
814
815       void Set( const vector< _Node* >& nodes, const vector< int >& quant = vector< int >() )
816       { _nodes.assign( nodes.begin(), nodes.end() ); _quantities = quant; }
817
818       void Set( _Node** nodes, int nb )
819       { _nodes.assign( nodes, nodes + nb ); }
820
821       void SetNext( _volumeDef* vd )
822       { if ( _next ) { _next->SetNext( vd ); } else { _next = vd; }}
823     };
824
825     // topology of a hexahedron
826     int   _nodeShift[8];
827     _Node _hexNodes [8];
828     _Link _hexLinks [12];
829     _Face _hexQuads [6];
830
831     // faces resulted from hexahedron intersection
832     vector< _Face > _polygons;
833
834     // intresections with EDGEs
835     vector< const E_IntersectPoint* > _eIntPoints;
836
837     // additional nodes created at intersection points
838     vector< _Node > _intNodes;
839
840     // nodes inside the hexahedron (at VERTEXes)
841     vector< _Node* > _vIntNodes;
842
843     // computed volume elements
844     _volumeDef _volumeDefs;
845
846     Grid*       _grid;
847     double      _sideLength[3];
848     int         _nbCornerNodes, _nbFaceIntNodes, _nbBndNodes;
849     int         _origNodeInd; // index of _hexNodes[0] node within the _grid
850     size_t      _i,_j,_k;
851     bool        _hasTooSmall;
852
853 #ifdef _DEBUG_
854     int         _cellID;
855 #endif
856
857   public:
858     Hexahedron(Grid* grid);
859     int MakeElements(SMESH_MesherHelper&                      helper,
860                      const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap);
861     void ComputeElements( const Solid* solid = 0, int solidIndex = -1 );
862
863   private:
864     Hexahedron(const Hexahedron& other, size_t i, size_t j, size_t k, int cellID );
865     void init( size_t i, size_t j, size_t k, const Solid* solid=0 );
866     void init( size_t i );
867     void setIJK( size_t i );
868     bool compute( const Solid* solid, const IsInternalFlag intFlag );
869     vector< TGeomID > getSolids();
870     bool isCutByInternalFace( IsInternalFlag & maxFlag );
871     void addEdges(SMESH_MesherHelper&                      helper,
872                   vector< Hexahedron* >&                   intersectedHex,
873                   const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap);
874     gp_Pnt findIntPoint( double u1, double proj1, double u2, double proj2,
875                          double proj, BRepAdaptor_Curve& curve,
876                          const gp_XYZ& axis, const gp_XYZ& origin );
877     int  getEntity( const E_IntersectPoint* ip, int* facets, int& sub );
878     bool addIntersection( const E_IntersectPoint* ip,
879                           vector< Hexahedron* >&  hexes,
880                           int ijk[], int dIJK[] );
881     bool findChain( _Node* n1, _Node* n2, _Face& quad, vector<_Node*>& chainNodes );
882     bool closePolygon( _Face* polygon, vector<_Node*>& chainNodes ) const;
883     bool findChainOnEdge( const vector< _OrientedLink >& splits,
884                           const _OrientedLink&           prevSplit,
885                           const _OrientedLink&           avoidSplit,
886                           size_t &                       iS,
887                           _Face&                         quad,
888                           vector<_Node*>&                chn);
889     int  addVolumes(SMESH_MesherHelper& helper );
890     void addFaces( SMESH_MesherHelper&                       helper,
891                    const vector< const SMDS_MeshElement* > & boundaryVolumes );
892     void addSegments( SMESH_MesherHelper&                      helper,
893                       const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap );
894     void getVolumes( vector< const SMDS_MeshElement* > & volumes );
895     void getBoundaryElems( vector< const SMDS_MeshElement* > & boundaryVolumes );
896     TGeomID getAnyFace() const;
897     void cutByExtendedInternal( std::vector< Hexahedron* >& hexes,
898                                 const TColStd_MapOfInteger& intEdgeIDs );
899     gp_Pnt mostDistantInternalPnt( int hexIndex, const gp_Pnt& p1, const gp_Pnt& p2 );
900     bool isOutPoint( _Link& link, int iP, SMESH_MesherHelper& helper, const Solid* solid ) const;
901     void sortVertexNodes(vector<_Node*>& nodes, _Node* curNode, TGeomID face);
902     bool isInHole() const;
903     bool hasStrangeEdge() const;
904     bool checkPolyhedronSize( bool isCutByInternalFace ) const;
905     bool addHexa ();
906     bool addTetra();
907     bool addPenta();
908     bool addPyra ();
909     bool debugDumpLink( _Link* link );
910     _Node* findEqualNode( vector< _Node* >&       nodes,
911                           const E_IntersectPoint* ip,
912                           const double            tol2 )
913     {
914       for ( size_t i = 0; i < nodes.size(); ++i )
915         if ( nodes[i]->EdgeIntPnt() == ip ||
916              nodes[i]->Point().SquareDistance( ip->_point ) <= tol2 )
917           return nodes[i];
918       return 0;
919     }
920     bool isImplementEdges() const { return _grid->_edgeIntPool.nbElements(); }
921     bool isOutParam(const double uvw[3]) const;
922
923     typedef boost::container::flat_map< TGeomID, size_t > TID2Nb;
924     static void insertAndIncrement( TGeomID id, TID2Nb& id2nbMap )
925     {
926       TID2Nb::value_type s0( id, 0 );
927       TID2Nb::iterator id2nb = id2nbMap.insert( s0 ).first;
928       id2nb->second++;
929     }
930   };
931
932 #ifdef WITH_TBB
933   // --------------------------------------------------------------------------
934   /*!
935    * \brief Hexahedron computing volumes in one thread
936    */
937   struct ParallelHexahedron
938   {
939     vector< Hexahedron* >& _hexVec;
940     ParallelHexahedron( vector< Hexahedron* >& hv ): _hexVec(hv) {}
941     void operator() ( const tbb::blocked_range<size_t>& r ) const
942     {
943       for ( size_t i = r.begin(); i != r.end(); ++i )
944         if ( Hexahedron* hex = _hexVec[ i ] )
945           hex->ComputeElements();
946     }
947   };
948   // --------------------------------------------------------------------------
949   /*!
950    * \brief Structure intersecting certain nb of faces with GridLine's in one thread
951    */
952   struct ParallelIntersector
953   {
954     vector< FaceGridIntersector >& _faceVec;
955     ParallelIntersector( vector< FaceGridIntersector >& faceVec): _faceVec(faceVec){}
956     void operator() ( const tbb::blocked_range<size_t>& r ) const
957     {
958       for ( size_t i = r.begin(); i != r.end(); ++i )
959         _faceVec[i].Intersect();
960     }
961   };
962 #endif
963
964   //=============================================================================
965   // Implementation of internal utils
966   //=============================================================================
967   /*!
968    * \brief adjust \a i to have \a val between values[i] and values[i+1]
969    */
970   inline void locateValue( int & i, double val, const vector<double>& values,
971                            int& di, double tol )
972   {
973     //val += values[0]; // input \a val is measured from 0.
974     if ( i > (int) values.size()-2 )
975       i = values.size()-2;
976     else
977       while ( i+2 < (int) values.size() && val > values[ i+1 ])
978         ++i;
979     while ( i > 0 && val < values[ i ])
980       --i;
981
982     if ( i > 0 && val - values[ i ] < tol )
983       di = -1;
984     else if ( i+2 < (int) values.size() && values[ i+1 ] - val < tol )
985       di = 1;
986     else
987       di = 0;
988   }
989   //=============================================================================
990   /*
991    * Remove coincident intersection points
992    */
993   void GridLine::RemoveExcessIntPoints( const double tol )
994   {
995     if ( _intPoints.size() < 2 ) return;
996
997     set< Transition > tranSet;
998     multiset< F_IntersectPoint >::iterator ip1, ip2 = _intPoints.begin();
999     while ( ip2 != _intPoints.end() )
1000     {
1001       tranSet.clear();
1002       ip1 = ip2++;
1003       while ( ip2 != _intPoints.end() && ip2->_paramOnLine - ip1->_paramOnLine <= tol )
1004       {
1005         tranSet.insert( ip1->_transition );
1006         tranSet.insert( ip2->_transition );
1007         ip2->Add( ip1->_faceIDs );
1008         _intPoints.erase( ip1 );
1009         ip1 = ip2++;
1010       }
1011       if ( tranSet.size() > 1 ) // points with different transition coincide
1012       {
1013         bool isIN  = tranSet.count( Trans_IN );
1014         bool isOUT = tranSet.count( Trans_OUT );
1015         if ( isIN && isOUT )
1016           (*ip1)._transition = Trans_TANGENT;
1017         else
1018           (*ip1)._transition = isIN ? Trans_IN : Trans_OUT;
1019       }
1020     }
1021   }
1022   //================================================================================
1023   /*
1024    * Return ID of SOLID for nodes before the given intersection point
1025    */
1026   TGeomID GridLine::GetSolidIDBefore( multiset< F_IntersectPoint >::iterator ip,
1027                                       const TGeomID                          prevID,
1028                                       const Geometry&                        geom )
1029   {
1030     if ( ip == _intPoints.begin() )
1031       return 0;
1032
1033     if ( geom.IsOneSolid() )
1034     {
1035       bool isOut = true;
1036       switch ( ip->_transition ) {
1037       case Trans_IN:      isOut = true;            break;
1038       case Trans_OUT:     isOut = false;           break;
1039       case Trans_TANGENT: isOut = ( prevID == 0 ); break;
1040       case Trans_APEX:
1041       {
1042         // singularity point (apex of a cone)
1043         multiset< F_IntersectPoint >::iterator ipBef = ip, ipAft = ++ip;
1044         if ( ipAft == _intPoints.end() )
1045           isOut = false;
1046         else
1047         {
1048           --ipBef;
1049           if ( ipBef->_transition != ipAft->_transition )
1050             isOut = ( ipBef->_transition == Trans_OUT );
1051           else
1052             isOut = ( ipBef->_transition != Trans_OUT );
1053         }
1054         break;
1055       }
1056       case Trans_INTERNAL: isOut = false;
1057       default:;
1058       }
1059       return isOut ? 0 : geom._soleSolid.ID();
1060     }
1061
1062     const vector< TGeomID >& solids = geom._solidIDsByShapeID[ ip->_faceIDs[ 0 ]];
1063
1064     --ip;
1065     if ( ip->_transition == Trans_INTERNAL )
1066       return prevID;
1067
1068     const vector< TGeomID >& solidsBef = geom._solidIDsByShapeID[ ip->_faceIDs[ 0 ]];
1069
1070     if ( ip->_transition == Trans_IN ||
1071          ip->_transition == Trans_OUT )
1072     {
1073       if ( solidsBef.size() == 1 )
1074         return ( solidsBef[0] == prevID ) ? 0 : solidsBef[0];
1075
1076       return solidsBef[ solidsBef[0] == prevID ];
1077     }
1078
1079     if ( solidsBef.size() == 1 )
1080       return solidsBef[0];
1081
1082     for ( size_t i = 0; i < solids.size(); ++i )
1083     {
1084       vector< TGeomID >::const_iterator it =
1085         std::find( solidsBef.begin(), solidsBef.end(), solids[i] );
1086       if ( it != solidsBef.end() )
1087         return solids[i];
1088     }
1089     return 0;
1090   }
1091   //================================================================================
1092   /*
1093    * Adds face IDs
1094    */
1095   void B_IntersectPoint::Add( const vector< TGeomID >& fIDs,
1096                               const SMDS_MeshNode*     n) const
1097   {
1098     if ( _faceIDs.empty() )
1099       _faceIDs = fIDs;
1100     else
1101       for ( size_t i = 0; i < fIDs.size(); ++i )
1102       {
1103         vector< TGeomID >::iterator it =
1104           std::find( _faceIDs.begin(), _faceIDs.end(), fIDs[i] );
1105         if ( it == _faceIDs.end() )
1106           _faceIDs.push_back( fIDs[i] );
1107       }
1108     if ( !_node )
1109       _node = n;
1110   }
1111   //================================================================================
1112   /*
1113    * Returns index of a common face if any, else zero
1114    */
1115   int B_IntersectPoint::HasCommonFace( const B_IntersectPoint * other, int avoidFace ) const
1116   {
1117     if ( other )
1118       for ( size_t i = 0; i < other->_faceIDs.size(); ++i )
1119         if ( avoidFace != other->_faceIDs[i] &&
1120              IsOnFace   ( other->_faceIDs[i] ))
1121           return other->_faceIDs[i];
1122     return 0;
1123   }
1124   //================================================================================
1125   /*
1126    * Returns \c true if \a faceID in in this->_faceIDs
1127    */
1128   bool B_IntersectPoint::IsOnFace( int faceID ) const // returns true if faceID is found
1129   {
1130     vector< TGeomID >::const_iterator it =
1131       std::find( _faceIDs.begin(), _faceIDs.end(), faceID );
1132     return ( it != _faceIDs.end() );
1133   }
1134   //================================================================================
1135   /*
1136    * OneOfSolids initialization
1137    */
1138   void OneOfSolids::Init( const TopoDS_Shape& solid,
1139                           TopAbs_ShapeEnum    subType,
1140                           const SMESHDS_Mesh* mesh )
1141   {
1142     SetID( mesh->ShapeToIndex( solid ));
1143
1144     if ( subType == TopAbs_FACE )
1145       SetHasInternalFaces( false );
1146
1147     for ( TopExp_Explorer sub( solid, subType ); sub.More(); sub.Next() )
1148     {
1149       _subIDs.Add( mesh->ShapeToIndex( sub.Current() ));
1150       if ( subType == TopAbs_FACE )
1151       {
1152         _faces.Add( sub.Current() );
1153         if ( sub.Current().Orientation() == TopAbs_INTERNAL )
1154           SetHasInternalFaces( true );
1155
1156         TGeomID faceID = mesh->ShapeToIndex( sub.Current() );
1157         if ( sub.Current().Orientation() == TopAbs_INTERNAL ||
1158              sub.Current().Orientation() == mesh->IndexToShape( faceID ).Orientation() )
1159           _outFaceIDs.Add( faceID );
1160       }
1161     }
1162   }
1163   //================================================================================
1164   /*
1165    * Return an iterator on GridLine's in a given direction
1166    */
1167   LineIndexer Grid::GetLineIndexer(size_t iDir) const
1168   {
1169     const size_t indices[] = { 1,2,0, 0,2,1, 0,1,2 };
1170     const string s      [] = { "X", "Y", "Z" };
1171     LineIndexer li( _coords[0].size(),  _coords[1].size(),    _coords[2].size(),
1172                     indices[iDir*3],    indices[iDir*3+1],    indices[iDir*3+2],
1173                     s[indices[iDir*3]], s[indices[iDir*3+1]], s[indices[iDir*3+2]]);
1174     return li;
1175   }
1176   //=============================================================================
1177   /*
1178    * Creates GridLine's of the grid
1179    */
1180   void Grid::SetCoordinates(const vector<double>& xCoords,
1181                             const vector<double>& yCoords,
1182                             const vector<double>& zCoords,
1183                             const double*         axesDirs,
1184                             const Bnd_Box&        shapeBox)
1185   {
1186     _coords[0] = xCoords;
1187     _coords[1] = yCoords;
1188     _coords[2] = zCoords;
1189
1190     _axes[0].SetCoord( axesDirs[0],
1191                        axesDirs[1],
1192                        axesDirs[2]);
1193     _axes[1].SetCoord( axesDirs[3],
1194                        axesDirs[4],
1195                        axesDirs[5]);
1196     _axes[2].SetCoord( axesDirs[6],
1197                        axesDirs[7],
1198                        axesDirs[8]);
1199     _axes[0].Normalize();
1200     _axes[1].Normalize();
1201     _axes[2].Normalize();
1202
1203     _invB.SetCols( _axes[0], _axes[1], _axes[2] );
1204     _invB.Invert();
1205
1206     // compute tolerance
1207     _minCellSize = Precision::Infinite();
1208     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
1209     {
1210       for ( size_t i = 1; i < _coords[ iDir ].size(); ++i )
1211       {
1212         double cellLen = _coords[ iDir ][ i ] - _coords[ iDir ][ i-1 ];
1213         if ( cellLen < _minCellSize )
1214           _minCellSize = cellLen;
1215       }
1216     }
1217     if ( _minCellSize < Precision::Confusion() )
1218       throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
1219                                 SMESH_Comment("Too small cell size: ") << _minCellSize );
1220     _tol = _minCellSize / 1000.;
1221
1222     // attune grid extremities to shape bounding box
1223
1224     double sP[6]; // aXmin, aYmin, aZmin, aXmax, aYmax, aZmax
1225     shapeBox.Get(sP[0],sP[1],sP[2],sP[3],sP[4],sP[5]);
1226     double* cP[6] = { &_coords[0].front(), &_coords[1].front(), &_coords[2].front(),
1227                       &_coords[0].back(),  &_coords[1].back(),  &_coords[2].back() };
1228     for ( int i = 0; i < 6; ++i )
1229       if ( fabs( sP[i] - *cP[i] ) < _tol )
1230         *cP[i] = sP[i];// + _tol/1000. * ( i < 3 ? +1 : -1 );
1231
1232     for ( int iDir = 0; iDir < 3; ++iDir )
1233     {
1234       if ( _coords[iDir][0] - sP[iDir] > _tol )
1235       {
1236         _minCellSize = Min( _minCellSize, _coords[iDir][0] - sP[iDir] );
1237         _coords[iDir].insert( _coords[iDir].begin(), sP[iDir] + _tol/1000.);
1238       }
1239       if ( sP[iDir+3] - _coords[iDir].back() > _tol  )
1240       {
1241         _minCellSize = Min( _minCellSize, sP[iDir+3] - _coords[iDir].back() );
1242         _coords[iDir].push_back( sP[iDir+3] - _tol/1000.);
1243       }
1244     }
1245     _tol = _minCellSize / 1000.;
1246
1247     _origin = ( _coords[0][0] * _axes[0] +
1248                 _coords[1][0] * _axes[1] +
1249                 _coords[2][0] * _axes[2] );
1250
1251     // create lines
1252     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
1253     {
1254       LineIndexer li = GetLineIndexer( iDir );
1255       _lines[iDir].resize( li.NbLines() );
1256       double len = _coords[ iDir ].back() - _coords[iDir].front();
1257       for ( ; li.More(); ++li )
1258       {
1259         GridLine& gl = _lines[iDir][ li.LineIndex() ];
1260         gl._line.SetLocation( _coords[0][li.I()] * _axes[0] +
1261                               _coords[1][li.J()] * _axes[1] +
1262                               _coords[2][li.K()] * _axes[2] );
1263         gl._line.SetDirection( _axes[ iDir ]);
1264         gl._length = len;
1265       }
1266     }
1267   }
1268   //================================================================================
1269   /*
1270    * Return local ID of shape
1271    */
1272   TGeomID Grid::ShapeID( const TopoDS_Shape& s ) const
1273   {
1274     return _helper->GetMeshDS()->ShapeToIndex( s );
1275   }
1276   //================================================================================
1277   /*
1278    * Return a shape by its local ID
1279    */
1280   const TopoDS_Shape& Grid::Shape( TGeomID id ) const
1281   {
1282     return _helper->GetMeshDS()->IndexToShape( id );
1283   }
1284   //================================================================================
1285   /*
1286    * Initialize _geometry
1287    */
1288   void Grid::InitGeometry( const TopoDS_Shape& theShapeToMesh )
1289   {
1290     SMESH_Mesh* mesh = _helper->GetMesh();
1291
1292     _geometry._mainShape = theShapeToMesh;
1293     _geometry._extIntFaceID = mesh->GetMeshDS()->MaxShapeIndex() * 100;
1294     _geometry._soleSolid.SetID( 0 );
1295     _geometry._soleSolid.SetHasInternalFaces( false );
1296
1297     InitClassifier( theShapeToMesh, TopAbs_VERTEX, _geometry._vertexClassifier );
1298     InitClassifier( theShapeToMesh, TopAbs_EDGE  , _geometry._edgeClassifier );
1299
1300     TopExp_Explorer solidExp( theShapeToMesh, TopAbs_SOLID );
1301
1302     bool isSeveralSolids = false;
1303     if ( _toConsiderInternalFaces ) // check nb SOLIDs
1304     {
1305       solidExp.Next();
1306       isSeveralSolids = solidExp.More();
1307       _toConsiderInternalFaces = isSeveralSolids;
1308       solidExp.ReInit();
1309
1310       if ( !isSeveralSolids ) // look for an internal FACE
1311       {
1312         TopExp_Explorer fExp( theShapeToMesh, TopAbs_FACE );
1313         for ( ; fExp.More() &&  !_toConsiderInternalFaces; fExp.Next() )
1314           _toConsiderInternalFaces = ( fExp.Current().Orientation() == TopAbs_INTERNAL );
1315
1316         _geometry._soleSolid.SetHasInternalFaces( _toConsiderInternalFaces );
1317         _geometry._soleSolid.SetID( ShapeID( solidExp.Current() ));
1318       }
1319       else // fill Geometry::_solidByID
1320       {
1321         for ( ; solidExp.More(); solidExp.Next() )
1322         {
1323           OneOfSolids & solid = _geometry._solidByID[ ShapeID( solidExp.Current() )];
1324           solid.Init( solidExp.Current(), TopAbs_FACE,   mesh->GetMeshDS() );
1325           solid.Init( solidExp.Current(), TopAbs_EDGE,   mesh->GetMeshDS() );
1326           solid.Init( solidExp.Current(), TopAbs_VERTEX, mesh->GetMeshDS() );
1327         }
1328       }
1329     }
1330     else
1331     {
1332       _geometry._soleSolid.SetID( ShapeID( solidExp.Current() ));
1333     }
1334
1335     if ( !_toCreateFaces )
1336     {
1337       int nbSolidsGlobal = _helper->Count( mesh->GetShapeToMesh(), TopAbs_SOLID, false );
1338       int nbSolidsLocal  = _helper->Count( theShapeToMesh,         TopAbs_SOLID, false );
1339       _toCreateFaces = ( nbSolidsLocal < nbSolidsGlobal );
1340     }
1341
1342     TopTools_IndexedMapOfShape faces;
1343     if ( _toCreateFaces || isSeveralSolids )
1344       TopExp::MapShapes( theShapeToMesh, TopAbs_FACE, faces );
1345
1346     // find boundary FACEs on boundary of mesh->ShapeToMesh()
1347     if ( _toCreateFaces )
1348       for ( int i = 1; i <= faces.Size(); ++i )
1349         if ( faces(i).Orientation() != TopAbs_INTERNAL &&
1350              _helper->NbAncestors( faces(i), *mesh, TopAbs_SOLID ) == 1 )
1351         {
1352           _geometry._boundaryFaces.Add( ShapeID( faces(i) ));
1353         }
1354
1355     if ( isSeveralSolids )
1356       for ( int i = 1; i <= faces.Size(); ++i )
1357       {
1358         SetSolidFather( faces(i), theShapeToMesh );
1359         for ( TopExp_Explorer eExp( faces(i), TopAbs_EDGE ); eExp.More(); eExp.Next() )
1360         {
1361           const TopoDS_Edge& edge = TopoDS::Edge( eExp.Current() );
1362           SetSolidFather( edge, theShapeToMesh );
1363           SetSolidFather( _helper->IthVertex( 0, edge ), theShapeToMesh );
1364           SetSolidFather( _helper->IthVertex( 1, edge ), theShapeToMesh );
1365         }
1366       }
1367     return;
1368   }
1369   //================================================================================
1370   /*
1371    * Store ID of SOLID as father of its child shape ID
1372    */
1373   void Grid::SetSolidFather( const TopoDS_Shape& s, const TopoDS_Shape& theShapeToMesh )
1374   {
1375     if ( _geometry._solidIDsByShapeID.empty() )
1376       _geometry._solidIDsByShapeID.resize( _helper->GetMeshDS()->MaxShapeIndex() + 1 );
1377
1378     vector< TGeomID > & solidIDs = _geometry._solidIDsByShapeID[ ShapeID( s )];
1379     if ( !solidIDs.empty() )
1380       return;
1381     solidIDs.reserve(2);
1382     PShapeIteratorPtr solidIt = _helper->GetAncestors( s,
1383                                                        *_helper->GetMesh(),
1384                                                        TopAbs_SOLID,
1385                                                        & theShapeToMesh );
1386     while ( const TopoDS_Shape* solid = solidIt->next() )
1387       solidIDs.push_back( ShapeID( *solid ));
1388   }
1389   //================================================================================
1390   /*
1391    * Return IDs of solids given sub-shape belongs to
1392    */
1393   const vector< TGeomID > & Grid::GetSolidIDs( TGeomID subShapeID ) const
1394   {
1395     return _geometry._solidIDsByShapeID[ subShapeID ];
1396   }
1397   //================================================================================
1398   /*
1399    * Check if a sub-shape belongs to several SOLIDs
1400    */
1401   bool Grid::IsShared( TGeomID shapeID ) const
1402   {
1403     return !_geometry.IsOneSolid() && ( _geometry._solidIDsByShapeID[ shapeID ].size() > 1 );
1404   }
1405   //================================================================================
1406   /*
1407    * Check if any of FACEs belongs to several SOLIDs
1408    */
1409   bool Grid::IsAnyShared( const std::vector< TGeomID >& faceIDs ) const
1410   {
1411     for ( size_t i = 0; i < faceIDs.size(); ++i )
1412       if ( IsShared( faceIDs[ i ]))
1413         return true;
1414     return false;
1415   }
1416   //================================================================================
1417   /*
1418    * Return Solid by ID
1419    */
1420   Solid* Grid::GetSolid( TGeomID solidID )
1421   {
1422     if ( !solidID || _geometry.IsOneSolid() || _geometry._solidByID.empty() )
1423       return & _geometry._soleSolid;
1424
1425     return & _geometry._solidByID[ solidID ];
1426   }
1427   //================================================================================
1428   /*
1429    * Return OneOfSolids by ID
1430    */
1431   Solid* Grid::GetOneOfSolids( TGeomID solidID )
1432   {
1433     map< TGeomID, OneOfSolids >::iterator is2s = _geometry._solidByID.find( solidID );
1434     if ( is2s != _geometry._solidByID.end() )
1435       return & is2s->second;
1436
1437     return & _geometry._soleSolid;
1438   }
1439   //================================================================================
1440   /*
1441    * Check if transition on given FACE is correct for a given SOLID
1442    */
1443   bool Grid::IsCorrectTransition( TGeomID faceID, const Solid* solid )
1444   {
1445     if ( _geometry.IsOneSolid() )
1446       return true;
1447
1448     const vector< TGeomID >& solidIDs = _geometry._solidIDsByShapeID[ faceID ];
1449     return solidIDs[0] == solid->ID();
1450   }
1451
1452   //================================================================================
1453   /*
1454    * Assign to geometry a node at FACE intersection
1455    */
1456   void Grid::SetOnShape( const SMDS_MeshNode* n, const F_IntersectPoint& ip, bool unset )
1457   {
1458     TopoDS_Shape s;
1459     SMESHDS_Mesh* mesh = _helper->GetMeshDS();
1460     if ( ip._faceIDs.size() == 1 )
1461     {
1462       mesh->SetNodeOnFace( n, ip._faceIDs[0], ip._u, ip._v );
1463     }
1464     else if ( _geometry._vertexClassifier.IsSatisfy( n, &s ))
1465     {
1466       if ( unset ) mesh->UnSetNodeOnShape( n );
1467       mesh->SetNodeOnVertex( n, TopoDS::Vertex( s ));
1468     }
1469     else if ( _geometry._edgeClassifier.IsSatisfy( n, &s ))
1470     {
1471       if ( unset ) mesh->UnSetNodeOnShape( n );
1472       mesh->SetNodeOnEdge( n, TopoDS::Edge( s ));
1473     }
1474     else if ( ip._faceIDs.size() > 0 )
1475     {
1476       mesh->SetNodeOnFace( n, ip._faceIDs[0], ip._u, ip._v );
1477     }
1478     else if ( !unset && _geometry.IsOneSolid() )
1479     {
1480       mesh->SetNodeInVolume( n, _geometry._soleSolid.ID() );
1481     }
1482   }
1483   //================================================================================
1484   /*
1485    * Initialize a classifier
1486    */
1487   void Grid::InitClassifier( const TopoDS_Shape&        mainShape,
1488                              TopAbs_ShapeEnum           shapeType,
1489                              Controls::ElementsOnShape& classifier )
1490   {
1491     TopTools_IndexedMapOfShape shapes;
1492     TopExp::MapShapes( mainShape, shapeType, shapes );
1493
1494     TopoDS_Compound compound; BRep_Builder builder;
1495     builder.MakeCompound( compound );
1496     for ( int i = 1; i <= shapes.Size(); ++i )
1497       builder.Add( compound, shapes(i) );
1498
1499     classifier.SetMesh( _helper->GetMeshDS() );
1500     //classifier.SetTolerance( _tol ); // _tol is not initialised
1501     classifier.SetShape( compound, SMDSAbs_Node );
1502   }
1503
1504   //================================================================================
1505   /*
1506    * Return EDGEs with FACEs to implement into the mesh
1507    */
1508   void Grid::GetEdgesToImplement( map< TGeomID, vector< TGeomID > > & edge2faceIDsMap,
1509                                   const TopoDS_Shape&                 shape,
1510                                   const vector< TopoDS_Shape >&       faces )
1511   {
1512     // check if there are strange EDGEs
1513     TopTools_IndexedMapOfShape faceMap;
1514     TopExp::MapShapes( _helper->GetMesh()->GetShapeToMesh(), TopAbs_FACE, faceMap );
1515     int nbFacesGlobal = faceMap.Size();
1516     faceMap.Clear( false );
1517     TopExp::MapShapes( shape, TopAbs_FACE, faceMap );
1518     int nbFacesLocal  = faceMap.Size();
1519     bool hasStrangeEdges = ( nbFacesGlobal > nbFacesLocal );
1520     if ( !_toAddEdges && !hasStrangeEdges )
1521       return; // no FACEs in contact with those meshed by other algo
1522
1523     for ( size_t i = 0; i < faces.size(); ++i )
1524     {
1525       _helper->SetSubShape( faces[i] );
1526       for ( TopExp_Explorer eExp( faces[i], TopAbs_EDGE ); eExp.More(); eExp.Next() )
1527       {
1528         const TopoDS_Edge& edge = TopoDS::Edge( eExp.Current() );
1529         if ( hasStrangeEdges )
1530         {
1531           bool hasStrangeFace = false;
1532           PShapeIteratorPtr faceIt = _helper->GetAncestors( edge, *_helper->GetMesh(), TopAbs_FACE);
1533           while ( const TopoDS_Shape* face = faceIt->next() )
1534             if (( hasStrangeFace = !faceMap.Contains( *face )))
1535               break;
1536           if ( !hasStrangeFace && !_toAddEdges )
1537             continue;
1538           _geometry._strangeEdges.Add( ShapeID( edge ));
1539           _geometry._strangeEdges.Add( ShapeID( _helper->IthVertex( 0, edge )));
1540           _geometry._strangeEdges.Add( ShapeID( _helper->IthVertex( 1, edge )));
1541         }
1542         if ( !SMESH_Algo::isDegenerated( edge ) &&
1543              !_helper->IsRealSeam( edge ))
1544         {
1545           edge2faceIDsMap[ ShapeID( edge )].push_back( ShapeID( faces[i] ));
1546         }
1547       }
1548     }
1549     return;
1550   }
1551
1552   //================================================================================
1553   /*
1554    * Computes coordinates of a point in the grid CS
1555    */
1556   void Grid::ComputeUVW(const gp_XYZ& P, double UVW[3])
1557   {
1558     gp_XYZ p = P * _invB;
1559     p.Coord( UVW[0], UVW[1], UVW[2] );
1560   }
1561   //================================================================================
1562   /*
1563    * Creates all nodes
1564    */
1565   void Grid::ComputeNodes(SMESH_MesherHelper& helper)
1566   {
1567     // state of each node of the grid relative to the geometry
1568     const size_t nbGridNodes = _coords[0].size() * _coords[1].size() * _coords[2].size();
1569     const TGeomID undefID = 1e+9;
1570     vector< TGeomID > shapeIDVec( nbGridNodes, undefID );
1571     _nodes.resize( nbGridNodes, 0 );
1572     _gridIntP.resize( nbGridNodes, NULL );
1573
1574     SMESHDS_Mesh* mesh = helper.GetMeshDS();
1575
1576     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
1577     {
1578       LineIndexer li = GetLineIndexer( iDir );
1579
1580       // find out a shift of node index while walking along a GridLine in this direction
1581       li.SetIndexOnLine( 0 );
1582       size_t nIndex0 = NodeIndex( li.I(), li.J(), li.K() );
1583       li.SetIndexOnLine( 1 );
1584       const size_t nShift = NodeIndex( li.I(), li.J(), li.K() ) - nIndex0;
1585       
1586       const vector<double> & coords = _coords[ iDir ];
1587       for ( ; li.More(); ++li ) // loop on lines in iDir
1588       {
1589         li.SetIndexOnLine( 0 );
1590         nIndex0 = NodeIndex( li.I(), li.J(), li.K() );
1591
1592         GridLine& line = _lines[ iDir ][ li.LineIndex() ];
1593         const gp_XYZ lineLoc = line._line.Location().XYZ();
1594         const gp_XYZ lineDir = line._line.Direction().XYZ();
1595
1596         line.RemoveExcessIntPoints( _tol );
1597         multiset< F_IntersectPoint >&     intPnts = line._intPoints;
1598         multiset< F_IntersectPoint >::iterator ip = intPnts.begin();
1599
1600         // Create mesh nodes at intersections with geometry
1601         // and set OUT state of nodes between intersections
1602
1603         TGeomID solidID = 0;
1604         const double* nodeCoord = & coords[0];
1605         const double* coord0    = nodeCoord;
1606         const double* coordEnd  = coord0 + coords.size();
1607         double nodeParam = 0;
1608         for ( ; ip != intPnts.end(); ++ip )
1609         {
1610           solidID = line.GetSolidIDBefore( ip, solidID, _geometry );
1611
1612           // set OUT state or just skip IN nodes before ip
1613           if ( nodeParam < ip->_paramOnLine - _tol )
1614           {
1615             while ( nodeParam < ip->_paramOnLine - _tol )
1616             {
1617               TGeomID & nodeShapeID = shapeIDVec[ nIndex0 + nShift * ( nodeCoord-coord0 ) ];
1618               nodeShapeID = Min( solidID, nodeShapeID );
1619               if ( ++nodeCoord <  coordEnd )
1620                 nodeParam = *nodeCoord - *coord0;
1621               else
1622                 break;
1623             }
1624             if ( nodeCoord == coordEnd ) break;
1625           }
1626           // create a mesh node on a GridLine at ip if it does not coincide with a grid node
1627           if ( nodeParam > ip->_paramOnLine + _tol )
1628           {
1629             gp_XYZ xyz = lineLoc + ip->_paramOnLine * lineDir;
1630             ip->_node = mesh->AddNode( xyz.X(), xyz.Y(), xyz.Z() );
1631             ip->_indexOnLine = nodeCoord-coord0-1;
1632             SetOnShape( ip->_node, *ip );
1633           }
1634           // create a mesh node at ip coincident with a grid node
1635           else
1636           {
1637             int nodeIndex = nIndex0 + nShift * ( nodeCoord-coord0 );
1638             if ( !_nodes[ nodeIndex ] )
1639             {
1640               gp_XYZ xyz = lineLoc + nodeParam * lineDir;
1641               _nodes   [ nodeIndex ] = mesh->AddNode( xyz.X(), xyz.Y(), xyz.Z() );
1642               //_gridIntP[ nodeIndex ] = & * ip;
1643               //SetOnShape( _nodes[ nodeIndex ], *ip );
1644             }
1645             if ( _gridIntP[ nodeIndex ] )
1646               _gridIntP[ nodeIndex ]->Add( ip->_faceIDs );
1647             else
1648               _gridIntP[ nodeIndex ] = & * ip;
1649             // ip->_node        = _nodes[ nodeIndex ]; -- to differ from ip on links
1650             ip->_indexOnLine = nodeCoord-coord0;
1651             if ( ++nodeCoord < coordEnd )
1652               nodeParam = *nodeCoord - *coord0;
1653           }
1654         }
1655         // set OUT state to nodes after the last ip
1656         for ( ; nodeCoord < coordEnd; ++nodeCoord )
1657           shapeIDVec[ nIndex0 + nShift * ( nodeCoord-coord0 ) ] = 0;
1658       }
1659     }
1660
1661     // Create mesh nodes at !OUT nodes of the grid
1662
1663     for ( size_t z = 0; z < _coords[2].size(); ++z )
1664       for ( size_t y = 0; y < _coords[1].size(); ++y )
1665         for ( size_t x = 0; x < _coords[0].size(); ++x )
1666         {
1667           size_t nodeIndex = NodeIndex( x, y, z );
1668           if ( !_nodes[ nodeIndex ] &&
1669                0 < shapeIDVec[ nodeIndex ] && shapeIDVec[ nodeIndex ] < undefID )
1670           {
1671             gp_XYZ xyz = ( _coords[0][x] * _axes[0] +
1672                            _coords[1][y] * _axes[1] +
1673                            _coords[2][z] * _axes[2] );
1674             _nodes[ nodeIndex ] = mesh->AddNode( xyz.X(), xyz.Y(), xyz.Z() );
1675             mesh->SetNodeInVolume( _nodes[ nodeIndex ], shapeIDVec[ nodeIndex ]);
1676           }
1677           else if ( _nodes[ nodeIndex ] && _gridIntP[ nodeIndex ] /*&&
1678                     !_nodes[ nodeIndex]->GetShapeID()*/ )
1679           {
1680             SetOnShape( _nodes[ nodeIndex ], *_gridIntP[ nodeIndex ]);
1681           }
1682         }
1683
1684 #ifdef _MY_DEBUG_
1685     // check validity of transitions
1686     const char* trName[] = { "TANGENT", "IN", "OUT", "APEX" };
1687     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
1688     {
1689       LineIndexer li = GetLineIndexer( iDir );
1690       for ( ; li.More(); ++li )
1691       {
1692         multiset< F_IntersectPoint >& intPnts = _lines[ iDir ][ li.LineIndex() ]._intPoints;
1693         if ( intPnts.empty() ) continue;
1694         if ( intPnts.size() == 1 )
1695         {
1696           if ( intPnts.begin()->_transition != Trans_TANGENT &&
1697                intPnts.begin()->_transition != Trans_APEX )
1698           throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
1699                                     SMESH_Comment("Wrong SOLE transition of GridLine (")
1700                                     << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
1701                                     << ") along " << li._nameConst
1702                                     << ": " << trName[ intPnts.begin()->_transition] );
1703         }
1704         else
1705         {
1706           if ( intPnts.begin()->_transition == Trans_OUT )
1707             throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
1708                                       SMESH_Comment("Wrong START transition of GridLine (")
1709                                       << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
1710                                       << ") along " << li._nameConst
1711                                       << ": " << trName[ intPnts.begin()->_transition ]);
1712           if ( intPnts.rbegin()->_transition == Trans_IN )
1713             throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
1714                                       SMESH_Comment("Wrong END transition of GridLine (")
1715                                       << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
1716                                       << ") along " << li._nameConst
1717                                     << ": " << trName[ intPnts.rbegin()->_transition ]);
1718         }
1719       }
1720     }
1721 #endif
1722   }
1723
1724   //=============================================================================
1725   /*
1726    * Intersects TopoDS_Face with all GridLine's
1727    */
1728   void FaceGridIntersector::Intersect()
1729   {
1730     FaceLineIntersector intersector;
1731     intersector._surfaceInt = GetCurveFaceIntersector();
1732     intersector._tol        = _grid->_tol;
1733     intersector._transOut   = _face.Orientation() == TopAbs_REVERSED ? Trans_IN : Trans_OUT;
1734     intersector._transIn    = _face.Orientation() == TopAbs_REVERSED ? Trans_OUT : Trans_IN;
1735
1736     typedef void (FaceLineIntersector::* PIntFun )(const GridLine& gridLine);
1737     PIntFun interFunction;
1738
1739     bool isDirect = true;
1740     BRepAdaptor_Surface surf( _face );
1741     switch ( surf.GetType() ) {
1742     case GeomAbs_Plane:
1743       intersector._plane = surf.Plane();
1744       interFunction = &FaceLineIntersector::IntersectWithPlane;
1745       isDirect = intersector._plane.Direct();
1746       break;
1747     case GeomAbs_Cylinder:
1748       intersector._cylinder = surf.Cylinder();
1749       interFunction = &FaceLineIntersector::IntersectWithCylinder;
1750       isDirect = intersector._cylinder.Direct();
1751       break;
1752     case GeomAbs_Cone:
1753       intersector._cone = surf.Cone();
1754       interFunction = &FaceLineIntersector::IntersectWithCone;
1755       //isDirect = intersector._cone.Direct();
1756       break;
1757     case GeomAbs_Sphere:
1758       intersector._sphere = surf.Sphere();
1759       interFunction = &FaceLineIntersector::IntersectWithSphere;
1760       isDirect = intersector._sphere.Direct();
1761       break;
1762     case GeomAbs_Torus:
1763       intersector._torus = surf.Torus();
1764       interFunction = &FaceLineIntersector::IntersectWithTorus;
1765       //isDirect = intersector._torus.Direct();
1766       break;
1767     default:
1768       interFunction = &FaceLineIntersector::IntersectWithSurface;
1769     }
1770     if ( !isDirect )
1771       std::swap( intersector._transOut, intersector._transIn );
1772
1773     _intersections.clear();
1774     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
1775     {
1776       if ( surf.GetType() == GeomAbs_Plane )
1777       {
1778         // check if all lines in this direction are parallel to a plane
1779         if ( intersector._plane.Axis().IsNormal( _grid->_lines[iDir][0]._line.Position(),
1780                                                  Precision::Angular()))
1781           continue;
1782         // find out a transition, that is the same for all lines of a direction
1783         gp_Dir plnNorm = intersector._plane.Axis().Direction();
1784         gp_Dir lineDir = _grid->_lines[iDir][0]._line.Direction();
1785         intersector._transition =
1786           ( plnNorm * lineDir < 0 ) ? intersector._transIn : intersector._transOut;
1787       }
1788       if ( surf.GetType() == GeomAbs_Cylinder )
1789       {
1790         // check if all lines in this direction are parallel to a cylinder
1791         if ( intersector._cylinder.Axis().IsParallel( _grid->_lines[iDir][0]._line.Position(),
1792                                                       Precision::Angular()))
1793           continue;
1794       }
1795
1796       // intersect the grid lines with the face
1797       for ( size_t iL = 0; iL < _grid->_lines[iDir].size(); ++iL )
1798       {
1799         GridLine& gridLine = _grid->_lines[iDir][iL];
1800         if ( _bndBox.IsOut( gridLine._line )) continue;
1801
1802         intersector._intPoints.clear();
1803         (intersector.*interFunction)( gridLine ); // <- intersection with gridLine
1804         for ( size_t i = 0; i < intersector._intPoints.size(); ++i )
1805           _intersections.push_back( make_pair( &gridLine, intersector._intPoints[i] ));
1806       }
1807     }
1808
1809     if ( _face.Orientation() == TopAbs_INTERNAL )
1810     {
1811       for ( size_t i = 0; i < _intersections.size(); ++i )
1812         if ( _intersections[i].second._transition == Trans_IN ||
1813              _intersections[i].second._transition == Trans_OUT )
1814         {
1815           _intersections[i].second._transition = Trans_INTERNAL;
1816         }
1817     }
1818     return;
1819   }
1820   //================================================================================
1821   /*
1822    * Return true if (_u,_v) is on the face
1823    */
1824   bool FaceLineIntersector::UVIsOnFace() const
1825   {
1826     TopAbs_State state = _surfaceInt->ClassifyUVPoint(gp_Pnt2d( _u,_v ));
1827     return ( state == TopAbs_IN || state == TopAbs_ON );
1828   }
1829   //================================================================================
1830   /*
1831    * Store an intersection if it is IN or ON the face
1832    */
1833   void FaceLineIntersector::addIntPoint(const bool toClassify)
1834   {
1835     if ( !toClassify || UVIsOnFace() )
1836     {
1837       F_IntersectPoint p;
1838       p._paramOnLine = _w;
1839       p._u           = _u;
1840       p._v           = _v;
1841       p._transition  = _transition;
1842       _intPoints.push_back( p );
1843     }
1844   }
1845   //================================================================================
1846   /*
1847    * Intersect a line with a plane
1848    */
1849   void FaceLineIntersector::IntersectWithPlane(const GridLine& gridLine)
1850   {
1851     IntAna_IntConicQuad linPlane( gridLine._line, _plane, Precision::Angular());
1852     _w = linPlane.ParamOnConic(1);
1853     if ( isParamOnLineOK( gridLine._length ))
1854     {
1855       ElSLib::Parameters(_plane, linPlane.Point(1) ,_u,_v);
1856       addIntPoint();
1857     }
1858   }
1859   //================================================================================
1860   /*
1861    * Intersect a line with a cylinder
1862    */
1863   void FaceLineIntersector::IntersectWithCylinder(const GridLine& gridLine)
1864   {
1865     IntAna_IntConicQuad linCylinder( gridLine._line, _cylinder );
1866     if ( linCylinder.IsDone() && linCylinder.NbPoints() > 0 )
1867     {
1868       _w = linCylinder.ParamOnConic(1);
1869       if ( linCylinder.NbPoints() == 1 )
1870         _transition = Trans_TANGENT;
1871       else
1872         _transition = _w < linCylinder.ParamOnConic(2) ? _transIn : _transOut;
1873       if ( isParamOnLineOK( gridLine._length ))
1874       {
1875         ElSLib::Parameters(_cylinder, linCylinder.Point(1) ,_u,_v);
1876         addIntPoint();
1877       }
1878       if ( linCylinder.NbPoints() > 1 )
1879       {
1880         _w = linCylinder.ParamOnConic(2);
1881         if ( isParamOnLineOK( gridLine._length ))
1882         {
1883           ElSLib::Parameters(_cylinder, linCylinder.Point(2) ,_u,_v);
1884           _transition = ( _transition == Trans_OUT ) ? Trans_IN : Trans_OUT;
1885           addIntPoint();
1886         }
1887       }
1888     }
1889   }
1890   //================================================================================
1891   /*
1892    * Intersect a line with a cone
1893    */
1894   void FaceLineIntersector::IntersectWithCone (const GridLine& gridLine)
1895   {
1896     IntAna_IntConicQuad linCone(gridLine._line,_cone);
1897     if ( !linCone.IsDone() ) return;
1898     gp_Pnt P;
1899     gp_Vec du, dv, norm;
1900     for ( int i = 1; i <= linCone.NbPoints(); ++i )
1901     {
1902       _w = linCone.ParamOnConic( i );
1903       if ( !isParamOnLineOK( gridLine._length )) continue;
1904       ElSLib::Parameters(_cone, linCone.Point(i) ,_u,_v);
1905       if ( UVIsOnFace() )
1906       {
1907         ElSLib::D1( _u, _v, _cone, P, du, dv );
1908         norm = du ^ dv;
1909         double normSize2 = norm.SquareMagnitude();
1910         if ( normSize2 > Precision::Angular() * Precision::Angular() )
1911         {
1912           double cos = norm.XYZ() * gridLine._line.Direction().XYZ();
1913           cos /= sqrt( normSize2 );
1914           if ( cos < -Precision::Angular() )
1915             _transition = _transIn;
1916           else if ( cos > Precision::Angular() )
1917             _transition = _transOut;
1918           else
1919             _transition = Trans_TANGENT;
1920         }
1921         else
1922         {
1923           _transition = Trans_APEX;
1924         }
1925         addIntPoint( /*toClassify=*/false);
1926       }
1927     }
1928   }
1929   //================================================================================
1930   /*
1931    * Intersect a line with a sphere
1932    */
1933   void FaceLineIntersector::IntersectWithSphere  (const GridLine& gridLine)
1934   {
1935     IntAna_IntConicQuad linSphere(gridLine._line,_sphere);
1936     if ( linSphere.IsDone() && linSphere.NbPoints() > 0 )
1937     {
1938       _w = linSphere.ParamOnConic(1);
1939       if ( linSphere.NbPoints() == 1 )
1940         _transition = Trans_TANGENT;
1941       else
1942         _transition = _w < linSphere.ParamOnConic(2) ? _transIn : _transOut;
1943       if ( isParamOnLineOK( gridLine._length ))
1944       {
1945         ElSLib::Parameters(_sphere, linSphere.Point(1) ,_u,_v);
1946         addIntPoint();
1947       }
1948       if ( linSphere.NbPoints() > 1 )
1949       {
1950         _w = linSphere.ParamOnConic(2);
1951         if ( isParamOnLineOK( gridLine._length ))
1952         {
1953           ElSLib::Parameters(_sphere, linSphere.Point(2) ,_u,_v);
1954           _transition = ( _transition == Trans_OUT ) ? Trans_IN : Trans_OUT;
1955           addIntPoint();
1956         }
1957       }
1958     }
1959   }
1960   //================================================================================
1961   /*
1962    * Intersect a line with a torus
1963    */
1964   void FaceLineIntersector::IntersectWithTorus   (const GridLine& gridLine)
1965   {
1966     IntAna_IntLinTorus linTorus(gridLine._line,_torus);
1967     if ( !linTorus.IsDone()) return;
1968     gp_Pnt P;
1969     gp_Vec du, dv, norm;
1970     for ( int i = 1; i <= linTorus.NbPoints(); ++i )
1971     {
1972       _w = linTorus.ParamOnLine( i );
1973       if ( !isParamOnLineOK( gridLine._length )) continue;
1974       linTorus.ParamOnTorus( i, _u,_v );
1975       if ( UVIsOnFace() )
1976       {
1977         ElSLib::D1( _u, _v, _torus, P, du, dv );
1978         norm = du ^ dv;
1979         double normSize = norm.Magnitude();
1980         double cos = norm.XYZ() * gridLine._line.Direction().XYZ();
1981         cos /= normSize;
1982         if ( cos < -Precision::Angular() )
1983           _transition = _transIn;
1984         else if ( cos > Precision::Angular() )
1985           _transition = _transOut;
1986         else
1987           _transition = Trans_TANGENT;
1988         addIntPoint( /*toClassify=*/false);
1989       }
1990     }
1991   }
1992   //================================================================================
1993   /*
1994    * Intersect a line with a non-analytical surface
1995    */
1996   void FaceLineIntersector::IntersectWithSurface (const GridLine& gridLine)
1997   {
1998     _surfaceInt->Perform( gridLine._line, 0.0, gridLine._length );
1999     if ( !_surfaceInt->IsDone() ) return;
2000     for ( int i = 1; i <= _surfaceInt->NbPnt(); ++i )
2001     {
2002       _transition = Transition( _surfaceInt->Transition( i ) );
2003       _w = _surfaceInt->WParameter( i );
2004       addIntPoint(/*toClassify=*/false);
2005     }
2006   }
2007   //================================================================================
2008   /*
2009    * check if its face can be safely intersected in a thread
2010    */
2011   bool FaceGridIntersector::IsThreadSafe(set< const Standard_Transient* >& noSafeTShapes) const
2012   {
2013     bool isSafe = true;
2014
2015     // check surface
2016     TopLoc_Location loc;
2017     Handle(Geom_Surface) surf = BRep_Tool::Surface( _face, loc );
2018     Handle(Geom_RectangularTrimmedSurface) ts =
2019       Handle(Geom_RectangularTrimmedSurface)::DownCast( surf );
2020     while( !ts.IsNull() ) {
2021       surf = ts->BasisSurface();
2022       ts = Handle(Geom_RectangularTrimmedSurface)::DownCast(surf);
2023     }
2024     if ( surf->IsKind( STANDARD_TYPE(Geom_BSplineSurface )) ||
2025          surf->IsKind( STANDARD_TYPE(Geom_BezierSurface )))
2026       if ( !noSafeTShapes.insert( _face.TShape().get() ).second )
2027         isSafe = false;
2028
2029     double f, l;
2030     TopExp_Explorer exp( _face, TopAbs_EDGE );
2031     for ( ; exp.More(); exp.Next() )
2032     {
2033       bool edgeIsSafe = true;
2034       const TopoDS_Edge& e = TopoDS::Edge( exp.Current() );
2035       // check 3d curve
2036       {
2037         Handle(Geom_Curve) c = BRep_Tool::Curve( e, loc, f, l);
2038         if ( !c.IsNull() )
2039         {
2040           Handle(Geom_TrimmedCurve) tc = Handle(Geom_TrimmedCurve)::DownCast(c);
2041           while( !tc.IsNull() ) {
2042             c = tc->BasisCurve();
2043             tc = Handle(Geom_TrimmedCurve)::DownCast(c);
2044           }
2045           if ( c->IsKind( STANDARD_TYPE(Geom_BSplineCurve )) ||
2046                c->IsKind( STANDARD_TYPE(Geom_BezierCurve )))
2047             edgeIsSafe = false;
2048         }
2049       }
2050       // check 2d curve
2051       if ( edgeIsSafe )
2052       {
2053         Handle(Geom2d_Curve) c2 = BRep_Tool::CurveOnSurface( e, surf, loc, f, l);
2054         if ( !c2.IsNull() )
2055         {
2056           Handle(Geom2d_TrimmedCurve) tc = Handle(Geom2d_TrimmedCurve)::DownCast(c2);
2057           while( !tc.IsNull() ) {
2058             c2 = tc->BasisCurve();
2059             tc = Handle(Geom2d_TrimmedCurve)::DownCast(c2);
2060           }
2061           if ( c2->IsKind( STANDARD_TYPE(Geom2d_BSplineCurve )) ||
2062                c2->IsKind( STANDARD_TYPE(Geom2d_BezierCurve )))
2063             edgeIsSafe = false;
2064         }
2065       }
2066       if ( !edgeIsSafe && !noSafeTShapes.insert( e.TShape().get() ).second )
2067         isSafe = false;
2068     }
2069     return isSafe;
2070   }
2071   //================================================================================
2072   /*!
2073    * \brief Creates topology of the hexahedron
2074    */
2075   Hexahedron::Hexahedron(Grid* grid)
2076     : _grid( grid ), _nbFaceIntNodes(0), _hasTooSmall( false )
2077   {
2078     _polygons.reserve(100); // to avoid reallocation;
2079
2080     //set nodes shift within grid->_nodes from the node 000 
2081     size_t dx = _grid->NodeIndexDX();
2082     size_t dy = _grid->NodeIndexDY();
2083     size_t dz = _grid->NodeIndexDZ();
2084     size_t i000 = 0;
2085     size_t i100 = i000 + dx;
2086     size_t i010 = i000 + dy;
2087     size_t i110 = i010 + dx;
2088     size_t i001 = i000 + dz;
2089     size_t i101 = i100 + dz;
2090     size_t i011 = i010 + dz;
2091     size_t i111 = i110 + dz;
2092     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V000 )] = i000;
2093     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V100 )] = i100;
2094     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V010 )] = i010;
2095     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V110 )] = i110;
2096     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V001 )] = i001;
2097     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V101 )] = i101;
2098     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V011 )] = i011;
2099     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V111 )] = i111;
2100
2101     vector< int > idVec;
2102     // set nodes to links
2103     for ( int linkID = SMESH_Block::ID_Ex00; linkID <= SMESH_Block::ID_E11z; ++linkID )
2104     {
2105       SMESH_Block::GetEdgeVertexIDs( linkID, idVec );
2106       _Link& link = _hexLinks[ SMESH_Block::ShapeIndex( linkID )];
2107       link._nodes[0] = &_hexNodes[ SMESH_Block::ShapeIndex( idVec[0] )];
2108       link._nodes[1] = &_hexNodes[ SMESH_Block::ShapeIndex( idVec[1] )];
2109     }
2110
2111     // set links to faces
2112     int interlace[4] = { 0, 3, 1, 2 }; // to walk by links around a face: { u0, 1v, u1, 0v }
2113     for ( int faceID = SMESH_Block::ID_Fxy0; faceID <= SMESH_Block::ID_F1yz; ++faceID )
2114     {
2115       SMESH_Block::GetFaceEdgesIDs( faceID, idVec );
2116       _Face& quad = _hexQuads[ SMESH_Block::ShapeIndex( faceID )];
2117       bool revFace = ( faceID == SMESH_Block::ID_Fxy0 ||
2118                        faceID == SMESH_Block::ID_Fx1z ||
2119                        faceID == SMESH_Block::ID_F0yz );
2120       quad._links.resize(4);
2121       vector<_OrientedLink>::iterator         frwLinkIt = quad._links.begin();
2122       vector<_OrientedLink>::reverse_iterator revLinkIt = quad._links.rbegin();
2123       for ( int i = 0; i < 4; ++i )
2124       {
2125         bool revLink = revFace;
2126         if ( i > 1 ) // reverse links u1 and v0
2127           revLink = !revLink;
2128         _OrientedLink& link = revFace ? *revLinkIt++ : *frwLinkIt++;
2129         link = _OrientedLink( & _hexLinks[ SMESH_Block::ShapeIndex( idVec[interlace[i]] )],
2130                               revLink );
2131       }
2132     }
2133   }
2134   //================================================================================
2135   /*!
2136    * \brief Copy constructor
2137    */
2138   Hexahedron::Hexahedron( const Hexahedron& other, size_t i, size_t j, size_t k, int cellID )
2139     :_grid( other._grid ), _nbFaceIntNodes(0), _i( i ), _j( j ), _k( k ), _hasTooSmall( false )
2140   {
2141     _polygons.reserve(100); // to avoid reallocation;
2142
2143     // copy topology
2144     for ( int i = 0; i < 8; ++i )
2145       _nodeShift[i] = other._nodeShift[i];
2146
2147     for ( int i = 0; i < 12; ++i )
2148     {
2149       const _Link& srcLink = other._hexLinks[ i ];
2150       _Link&       tgtLink = this->_hexLinks[ i ];
2151       tgtLink._nodes[0] = _hexNodes + ( srcLink._nodes[0] - other._hexNodes );
2152       tgtLink._nodes[1] = _hexNodes + ( srcLink._nodes[1] - other._hexNodes );
2153     }
2154
2155     for ( int i = 0; i < 6; ++i )
2156     {
2157       const _Face& srcQuad = other._hexQuads[ i ];
2158       _Face&       tgtQuad = this->_hexQuads[ i ];
2159       tgtQuad._links.resize(4);
2160       for ( int j = 0; j < 4; ++j )
2161       {
2162         const _OrientedLink& srcLink = srcQuad._links[ j ];
2163         _OrientedLink&       tgtLink = tgtQuad._links[ j ];
2164         tgtLink._reverse = srcLink._reverse;
2165         tgtLink._link    = _hexLinks + ( srcLink._link - other._hexLinks );
2166       }
2167     }
2168 #ifdef _DEBUG_
2169     _cellID = cellID;
2170 #endif
2171   }
2172
2173   //================================================================================
2174   /*!
2175    * \brief Return IDs of SOLIDs interfering with this Hexahedron
2176    */
2177   vector< TGeomID > Hexahedron::getSolids()
2178   {
2179     // count intersection points belonging to each SOLID
2180     TID2Nb id2NbPoints;
2181     id2NbPoints.reserve( 3 );
2182
2183     _origNodeInd = _grid->NodeIndex( _i,_j,_k );
2184     for ( int iN = 0; iN < 8; ++iN )
2185     {
2186       _hexNodes[iN]._node     = _grid->_nodes   [ _origNodeInd + _nodeShift[iN] ];
2187       _hexNodes[iN]._intPoint = _grid->_gridIntP[ _origNodeInd + _nodeShift[iN] ];
2188
2189       if ( _hexNodes[iN]._intPoint ) // intersection with a FACE
2190       {
2191         for ( size_t iF = 0; iF < _hexNodes[iN]._intPoint->_faceIDs.size(); ++iF )
2192         {
2193           const vector< TGeomID > & solidIDs =
2194             _grid->GetSolidIDs( _hexNodes[iN]._intPoint->_faceIDs[iF] );
2195           for ( size_t i = 0; i < solidIDs.size(); ++i )
2196             insertAndIncrement( solidIDs[i], id2NbPoints );
2197         }
2198       }
2199       else if ( _hexNodes[iN]._node ) // node inside a SOLID
2200       {
2201         insertAndIncrement( _hexNodes[iN]._node->GetShapeID(), id2NbPoints );
2202       }
2203     }
2204
2205     for ( int iL = 0; iL < 12; ++iL )
2206     {
2207       const _Link& link = _hexLinks[ iL ];
2208       for ( size_t iP = 0; iP < link._fIntPoints.size(); ++iP )
2209       {
2210         for ( size_t iF = 0; iF < link._fIntPoints[iP]->_faceIDs.size(); ++iF )
2211         {
2212           const vector< TGeomID > & solidIDs =
2213             _grid->GetSolidIDs( link._fIntPoints[iP]->_faceIDs[iF] );
2214           for ( size_t i = 0; i < solidIDs.size(); ++i )
2215             insertAndIncrement( solidIDs[i], id2NbPoints );
2216         }
2217       }
2218     }
2219
2220     for ( size_t iP = 0; iP < _eIntPoints.size(); ++iP )
2221     {
2222       const vector< TGeomID > & solidIDs = _grid->GetSolidIDs( _eIntPoints[iP]->_shapeID );
2223       for ( size_t i = 0; i < solidIDs.size(); ++i )
2224         insertAndIncrement( solidIDs[i], id2NbPoints );
2225     }
2226
2227     vector< TGeomID > solids; solids.reserve( id2NbPoints.size() );
2228     for ( TID2Nb::iterator id2nb = id2NbPoints.begin(); id2nb != id2NbPoints.end(); ++id2nb )
2229       if ( id2nb->second >= 3 )
2230         solids.push_back( id2nb->first );
2231
2232     return solids;
2233   }
2234
2235   //================================================================================
2236   /*!
2237    * \brief Count cuts by INTERNAL FACEs and set _Node::_isInternalFlags
2238    */
2239   bool Hexahedron::isCutByInternalFace( IsInternalFlag & maxFlag )
2240   {
2241     TID2Nb id2NbPoints;
2242     id2NbPoints.reserve( 3 );
2243
2244     for ( size_t iN = 0; iN < _intNodes.size(); ++iN )
2245       for ( size_t iF = 0; iF < _intNodes[iN]._intPoint->_faceIDs.size(); ++iF )
2246       {
2247         if ( _grid->IsInternal( _intNodes[iN]._intPoint->_faceIDs[iF]))
2248           insertAndIncrement( _intNodes[iN]._intPoint->_faceIDs[iF], id2NbPoints );
2249       }
2250     for ( size_t iN = 0; iN < 8; ++iN )
2251       if ( _hexNodes[iN]._intPoint )
2252         for ( size_t iF = 0; iF < _hexNodes[iN]._intPoint->_faceIDs.size(); ++iF )
2253         {
2254           if ( _grid->IsInternal( _hexNodes[iN]._intPoint->_faceIDs[iF]))
2255             insertAndIncrement( _hexNodes[iN]._intPoint->_faceIDs[iF], id2NbPoints );
2256         }
2257
2258     maxFlag = IS_NOT_INTERNAL;
2259     for ( TID2Nb::iterator id2nb = id2NbPoints.begin(); id2nb != id2NbPoints.end(); ++id2nb )
2260     {
2261       TGeomID        intFace = id2nb->first;
2262       IsInternalFlag intFlag = ( id2nb->second >= 3 ? IS_CUT_BY_INTERNAL_FACE : IS_INTERNAL );
2263       if ( intFlag > maxFlag )
2264         maxFlag = intFlag;
2265
2266       for ( size_t iN = 0; iN < _intNodes.size(); ++iN )
2267         if ( _intNodes[iN].IsOnFace( intFace ))
2268           _intNodes[iN].SetInternal( intFlag );
2269
2270       for ( size_t iN = 0; iN < 8; ++iN )
2271         if ( _hexNodes[iN].IsOnFace( intFace ))
2272           _hexNodes[iN].SetInternal( intFlag );
2273     }
2274
2275     return maxFlag;
2276   }
2277
2278   //================================================================================
2279   /*!
2280    * \brief Return any FACE interfering with this Hexahedron
2281    */
2282   TGeomID Hexahedron::getAnyFace() const
2283   {
2284     TID2Nb id2NbPoints;
2285     id2NbPoints.reserve( 3 );
2286
2287     for ( size_t iN = 0; iN < _intNodes.size(); ++iN )
2288       for ( size_t iF = 0; iF < _intNodes[iN]._intPoint->_faceIDs.size(); ++iF )
2289         insertAndIncrement( _intNodes[iN]._intPoint->_faceIDs[iF], id2NbPoints );
2290
2291     for ( size_t iN = 0; iN < 8; ++iN )
2292       if ( _hexNodes[iN]._intPoint )
2293         for ( size_t iF = 0; iF < _hexNodes[iN]._intPoint->_faceIDs.size(); ++iF )
2294           insertAndIncrement( _hexNodes[iN]._intPoint->_faceIDs[iF], id2NbPoints );
2295
2296     for ( unsigned int minNb = 3; minNb > 0; --minNb )
2297       for ( TID2Nb::iterator id2nb = id2NbPoints.begin(); id2nb != id2NbPoints.end(); ++id2nb )
2298         if ( id2nb->second >= minNb )
2299           return id2nb->first;
2300
2301     return 0;
2302   }
2303
2304   //================================================================================
2305   /*!
2306    * \brief Initializes IJK by Hexahedron index
2307    */
2308   void Hexahedron::setIJK( size_t iCell )
2309   {
2310     size_t iNbCell = _grid->_coords[0].size() - 1;
2311     size_t jNbCell = _grid->_coords[1].size() - 1;
2312     _i = iCell % iNbCell;
2313     _j = ( iCell % ( iNbCell * jNbCell )) / iNbCell;
2314     _k = iCell / iNbCell / jNbCell;
2315   }
2316
2317   //================================================================================
2318   /*!
2319    * \brief Initializes its data by given grid cell (countered from zero)
2320    */
2321   void Hexahedron::init( size_t iCell )
2322   {
2323     setIJK( iCell );
2324     init( _i, _j, _k );
2325   }
2326
2327   //================================================================================
2328   /*!
2329    * \brief Initializes its data by given grid cell nodes and intersections
2330    */
2331   void Hexahedron::init( size_t i, size_t j, size_t k, const Solid* solid )
2332   {
2333     _i = i; _j = j; _k = k;
2334
2335     if ( !solid )
2336       solid = _grid->GetSolid();
2337
2338     // set nodes of grid to nodes of the hexahedron and
2339     // count nodes at hexahedron corners located IN and ON geometry
2340     _nbCornerNodes = _nbBndNodes = 0;
2341     _origNodeInd   = _grid->NodeIndex( i,j,k );
2342     for ( int iN = 0; iN < 8; ++iN )
2343     {
2344       _hexNodes[iN]._isInternalFlags = 0;
2345
2346       _hexNodes[iN]._node     = _grid->_nodes   [ _origNodeInd + _nodeShift[iN] ];
2347       _hexNodes[iN]._intPoint = _grid->_gridIntP[ _origNodeInd + _nodeShift[iN] ];
2348
2349       if ( _hexNodes[iN]._node && !solid->Contains( _hexNodes[iN]._node->GetShapeID() ))
2350         _hexNodes[iN]._node = 0;
2351       if ( _hexNodes[iN]._intPoint && !solid->ContainsAny( _hexNodes[iN]._intPoint->_faceIDs ))
2352         _hexNodes[iN]._intPoint = 0;
2353
2354       _nbCornerNodes += bool( _hexNodes[iN]._node );
2355       _nbBndNodes    += bool( _hexNodes[iN]._intPoint );
2356     }
2357     _sideLength[0] = _grid->_coords[0][i+1] - _grid->_coords[0][i];
2358     _sideLength[1] = _grid->_coords[1][j+1] - _grid->_coords[1][j];
2359     _sideLength[2] = _grid->_coords[2][k+1] - _grid->_coords[2][k];
2360
2361     _intNodes.clear();
2362     _vIntNodes.clear();
2363
2364     if ( _nbFaceIntNodes + _eIntPoints.size()                  > 0 &&
2365          _nbFaceIntNodes + _eIntPoints.size() + _nbCornerNodes > 3)
2366     {
2367       _intNodes.reserve( 3 * _nbBndNodes + _nbFaceIntNodes + _eIntPoints.size() );
2368
2369       // this method can be called in parallel, so use own helper
2370       SMESH_MesherHelper helper( *_grid->_helper->GetMesh() );
2371
2372       // Create sub-links (_Link::_splits) by splitting links with _Link::_fIntPoints
2373       // ---------------------------------------------------------------
2374       _Link split;
2375       for ( int iLink = 0; iLink < 12; ++iLink )
2376       {
2377         _Link& link = _hexLinks[ iLink ];
2378         link._fIntNodes.clear();
2379         link._fIntNodes.reserve( link._fIntPoints.size() );
2380         for ( size_t i = 0; i < link._fIntPoints.size(); ++i )
2381           if ( solid->ContainsAny( link._fIntPoints[i]->_faceIDs ))
2382           {
2383             _intNodes.push_back( _Node( 0, link._fIntPoints[i] ));
2384             link._fIntNodes.push_back( & _intNodes.back() );
2385           }
2386
2387         link._splits.clear();
2388         split._nodes[ 0 ] = link._nodes[0];
2389         bool isOut = ( ! link._nodes[0]->Node() );
2390         bool checkTransition;
2391         for ( size_t i = 0; i < link._fIntNodes.size(); ++i )
2392         {
2393           const bool isGridNode = ( ! link._fIntNodes[i]->Node() );
2394           if ( !isGridNode ) // intersection non-coincident with a grid node
2395           {
2396             if ( split._nodes[ 0 ]->Node() && !isOut )
2397             {
2398               split._nodes[ 1 ] = link._fIntNodes[i];
2399               link._splits.push_back( split );
2400             }
2401             split._nodes[ 0 ] = link._fIntNodes[i];
2402             checkTransition = true;
2403           }
2404           else // FACE intersection coincident with a grid node (at link ends)
2405           {
2406             checkTransition = ( i == 0 && link._nodes[0]->Node() );
2407           }
2408           if ( checkTransition )
2409           {
2410             const vector< TGeomID >& faceIDs = link._fIntNodes[i]->_intPoint->_faceIDs;
2411             if ( _grid->IsInternal( faceIDs.back() ))
2412               isOut = false;
2413             else if ( faceIDs.size() > 1 || _eIntPoints.size() > 0 )
2414               isOut = isOutPoint( link, i, helper, solid );
2415             else
2416             {
2417               bool okTransi = _grid->IsCorrectTransition( faceIDs[0], solid );
2418               switch ( link._fIntNodes[i]->FaceIntPnt()->_transition ) {
2419               case Trans_OUT: isOut = okTransi;  break;
2420               case Trans_IN : isOut = !okTransi; break;
2421               default:
2422                 isOut = isOutPoint( link, i, helper, solid );
2423               }
2424             }
2425           }
2426         }
2427         if ( link._nodes[ 1 ]->Node() && split._nodes[ 0 ]->Node() && !isOut )
2428         {
2429           split._nodes[ 1 ] = link._nodes[1];
2430           link._splits.push_back( split );
2431         }
2432       }
2433
2434       // Create _Node's at intersections with EDGEs.
2435       // --------------------------------------------
2436       // 1) add this->_eIntPoints to _Face::_eIntNodes
2437       // 2) fill _intNodes and _vIntNodes
2438       //
2439       const double tol2 = _grid->_tol * _grid->_tol;
2440       int facets[3], nbFacets, subEntity;
2441
2442       for ( int iF = 0; iF < 6; ++iF )
2443         _hexQuads[ iF ]._eIntNodes.clear();
2444
2445       for ( size_t iP = 0; iP < _eIntPoints.size(); ++iP )
2446       {
2447         if ( !solid->ContainsAny( _eIntPoints[iP]->_faceIDs ))
2448           continue;
2449         nbFacets = getEntity( _eIntPoints[iP], facets, subEntity );
2450         _Node* equalNode = 0;
2451         switch( nbFacets ) {
2452         case 1: // in a _Face
2453         {
2454           _Face& quad = _hexQuads[ facets[0] - SMESH_Block::ID_FirstF ];
2455           equalNode = findEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
2456           if ( equalNode ) {
2457             equalNode->Add( _eIntPoints[ iP ] );
2458           }
2459           else {
2460             _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
2461             quad._eIntNodes.push_back( & _intNodes.back() );
2462           }
2463           break;
2464         }
2465         case 2: // on a _Link
2466         {
2467           _Link& link = _hexLinks[ subEntity - SMESH_Block::ID_FirstE ];
2468           if ( link._splits.size() > 0 )
2469           {
2470             equalNode = findEqualNode( link._fIntNodes, _eIntPoints[ iP ], tol2 );
2471             if ( equalNode )
2472               equalNode->Add( _eIntPoints[ iP ] );
2473             else if ( link._splits.size() == 1 &&
2474                       link._splits[0]._nodes[0] &&
2475                       link._splits[0]._nodes[1] )
2476               link._splits.clear(); // hex edge is divided by _eIntPoints[iP]
2477           }
2478           //else
2479           if ( !equalNode )
2480           {
2481             _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
2482             bool newNodeUsed = false;
2483             for ( int iF = 0; iF < 2; ++iF )
2484             {
2485               _Face& quad = _hexQuads[ facets[iF] - SMESH_Block::ID_FirstF ];
2486               equalNode = findEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
2487               if ( equalNode ) {
2488                 equalNode->Add( _eIntPoints[ iP ] );
2489               }
2490               else {
2491                 quad._eIntNodes.push_back( & _intNodes.back() );
2492                 newNodeUsed = true;
2493               }
2494             }
2495             if ( !newNodeUsed )
2496               _intNodes.pop_back();
2497           }
2498           break;
2499         }
2500         case 3: // at a corner
2501         {
2502           _Node& node = _hexNodes[ subEntity - SMESH_Block::ID_FirstV ];
2503           if ( node.Node() > 0 )
2504           {
2505             if ( node._intPoint )
2506               node._intPoint->Add( _eIntPoints[ iP ]->_faceIDs, _eIntPoints[ iP ]->_node );
2507           }
2508           else
2509           {
2510             _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
2511             for ( int iF = 0; iF < 3; ++iF )
2512             {
2513               _Face& quad = _hexQuads[ facets[iF] - SMESH_Block::ID_FirstF ];
2514               equalNode = findEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
2515               if ( equalNode ) {
2516                 equalNode->Add( _eIntPoints[ iP ] );
2517               }
2518               else {
2519                 quad._eIntNodes.push_back( & _intNodes.back() );
2520               }
2521             }
2522           }
2523           break;
2524         }
2525         } // switch( nbFacets )
2526
2527         if ( nbFacets == 0 ||
2528              _grid->ShapeType( _eIntPoints[ iP ]->_shapeID ) == TopAbs_VERTEX )
2529         {
2530           equalNode = findEqualNode( _vIntNodes, _eIntPoints[ iP ], tol2 );
2531           if ( equalNode ) {
2532             equalNode->Add( _eIntPoints[ iP ] );
2533           }
2534           else if ( nbFacets == 0 ) {
2535             if ( _intNodes.empty() || _intNodes.back().EdgeIntPnt() != _eIntPoints[ iP ])
2536               _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
2537             _vIntNodes.push_back( & _intNodes.back() );
2538           }
2539         }
2540       } // loop on _eIntPoints
2541     }
2542
2543     else if ( 3 < _nbCornerNodes && _nbCornerNodes < 8 ) // _nbFaceIntNodes == 0
2544     {
2545       _Link split;
2546       // create sub-links (_splits) of whole links
2547       for ( int iLink = 0; iLink < 12; ++iLink )
2548       {
2549         _Link& link = _hexLinks[ iLink ];
2550         link._splits.clear();
2551         if ( link._nodes[ 0 ]->Node() && link._nodes[ 1 ]->Node() )
2552         {
2553           split._nodes[ 0 ] = link._nodes[0];
2554           split._nodes[ 1 ] = link._nodes[1];
2555           link._splits.push_back( split );
2556         }
2557       }
2558     }
2559     return;
2560
2561   } // init( _i, _j, _k )
2562
2563   //================================================================================
2564   /*!
2565    * \brief Compute mesh volumes resulted from intersection of the Hexahedron
2566    */
2567   void Hexahedron::ComputeElements( const Solid* solid, int solidIndex )
2568   {
2569     if ( !solid )
2570     {
2571       solid = _grid->GetSolid();
2572       if ( !_grid->_geometry.IsOneSolid() )
2573       {
2574         vector< TGeomID > solidIDs = getSolids();
2575         if ( solidIDs.size() > 1 )
2576         {
2577           for ( size_t i = 0; i < solidIDs.size(); ++i )
2578           {
2579             solid = _grid->GetSolid( solidIDs[i] );
2580             ComputeElements( solid, i );
2581             if ( !_volumeDefs._nodes.empty() && i < solidIDs.size() - 1 )
2582               _volumeDefs.SetNext( new _volumeDef( _volumeDefs ));
2583           }
2584           return;
2585         }
2586         solid = _grid->GetSolid( solidIDs[0] );
2587       }
2588     }
2589
2590     init( _i, _j, _k, solid ); // get nodes and intersections from grid nodes and split links
2591
2592     int nbIntersections = _nbFaceIntNodes + _eIntPoints.size();
2593     if ( _nbCornerNodes + nbIntersections < 4 )
2594       return;
2595
2596     if ( _nbBndNodes == _nbCornerNodes && nbIntersections == 0 && isInHole() )
2597       return; // cell is in a hole
2598
2599     IsInternalFlag intFlag = IS_NOT_INTERNAL;
2600     if ( solid->HasInternalFaces() && this->isCutByInternalFace( intFlag ))
2601     {
2602       for ( _SplitIterator it( _hexLinks ); it.More(); it.Next() )
2603       {
2604         if ( compute( solid, intFlag ))
2605           _volumeDefs.SetNext( new _volumeDef( _volumeDefs ));
2606       }
2607     }
2608     else
2609     {
2610       if ( solidIndex >= 0 )
2611         intFlag = IS_CUT_BY_INTERNAL_FACE;
2612
2613       compute( solid, intFlag );
2614     }
2615   }
2616
2617   //================================================================================
2618   /*!
2619    * \brief Compute mesh volumes resulted from intersection of the Hexahedron
2620    */
2621   bool Hexahedron::compute( const Solid* solid, const IsInternalFlag intFlag )
2622   {
2623     _polygons.clear();
2624     _polygons.reserve( 20 );
2625
2626     for ( int iN = 0; iN < 8; ++iN )
2627       _hexNodes[iN]._usedInFace = 0;
2628
2629     // Create polygons from quadrangles
2630     // --------------------------------
2631
2632     vector< _OrientedLink > splits;
2633     vector<_Node*>          chainNodes;
2634     _Face*                  coplanarPolyg;
2635
2636     bool hasEdgeIntersections = !_eIntPoints.empty();
2637
2638     for ( int iF = 0; iF < 6; ++iF ) // loop on 6 sides of a hexahedron
2639     {
2640       _Face& quad = _hexQuads[ iF ] ;
2641
2642       _polygons.resize( _polygons.size() + 1 );
2643       _Face* polygon = &_polygons.back();
2644       polygon->_polyLinks.reserve( 20 );
2645
2646       splits.clear();
2647       for ( int iE = 0; iE < 4; ++iE ) // loop on 4 sides of a quadrangle
2648         for ( int iS = 0; iS < quad._links[ iE ].NbResultLinks(); ++iS )
2649           splits.push_back( quad._links[ iE ].ResultLink( iS ));
2650
2651       // add splits of links to a polygon and add _polyLinks to make
2652       // polygon's boundary closed
2653
2654       int nbSplits = splits.size();
2655       if (( nbSplits == 1 ) &&
2656           ( quad._eIntNodes.empty() ||
2657             splits[0].FirstNode()->IsLinked( splits[0].LastNode()->_intPoint )))
2658         //( quad._eIntNodes.empty() || _nbCornerNodes + nbIntersections > 6 ))
2659         nbSplits = 0;
2660
2661       for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
2662         if ( quad._eIntNodes[ iP ]->IsUsedInFace( polygon ))
2663           quad._eIntNodes[ iP ]->_usedInFace = 0;
2664
2665       size_t nbUsedEdgeNodes = 0;
2666       _Face* prevPolyg = 0; // polygon previously created from this quad
2667
2668       while ( nbSplits > 0 )
2669       {
2670         size_t iS = 0;
2671         while ( !splits[ iS ] )
2672           ++iS;
2673
2674         if ( !polygon->_links.empty() )
2675         {
2676           _polygons.resize( _polygons.size() + 1 );
2677           polygon = &_polygons.back();
2678           polygon->_polyLinks.reserve( 20 );
2679         }
2680         polygon->_links.push_back( splits[ iS ] );
2681         splits[ iS++ ]._link = 0;
2682         --nbSplits;
2683
2684         _Node* nFirst = polygon->_links.back().FirstNode();
2685         _Node *n1,*n2 = polygon->_links.back().LastNode();
2686         for ( ; nFirst != n2 && iS < splits.size(); ++iS )
2687         {
2688           _OrientedLink& split = splits[ iS ];
2689           if ( !split ) continue;
2690
2691           n1 = split.FirstNode();
2692           if ( n1 == n2 &&
2693                n1->_intPoint &&
2694                (( n1->_intPoint->_faceIDs.size() > 1 && isImplementEdges() ) ||
2695                 ( n1->_isInternalFlags )))
2696           {
2697             // n1 is at intersection with EDGE
2698             if ( findChainOnEdge( splits, polygon->_links.back(), split, iS, quad, chainNodes ))
2699             {
2700               for ( size_t i = 1; i < chainNodes.size(); ++i )
2701                 polygon->AddPolyLink( chainNodes[i-1], chainNodes[i], prevPolyg );
2702               if ( chainNodes.back() != n1 ) // not a partial cut by INTERNAL FACE
2703               {
2704                 prevPolyg = polygon;
2705                 n2 = chainNodes.back();
2706                 continue;
2707               }
2708             }
2709           }
2710           else if ( n1 != n2 )
2711           {
2712             // try to connect to intersections with EDGEs
2713             if ( quad._eIntNodes.size() > nbUsedEdgeNodes  &&
2714                  findChain( n2, n1, quad, chainNodes ))
2715             {
2716               for ( size_t i = 1; i < chainNodes.size(); ++i )
2717               {
2718                 polygon->AddPolyLink( chainNodes[i-1], chainNodes[i] );
2719                 nbUsedEdgeNodes += ( chainNodes[i]->IsUsedInFace( polygon ));
2720               }
2721               if ( chainNodes.back() != n1 )
2722               {
2723                 n2 = chainNodes.back();
2724                 --iS;
2725                 continue;
2726               }
2727             }
2728             // try to connect to a split ending on the same FACE
2729             else
2730             {
2731               _OrientedLink foundSplit;
2732               for ( size_t i = iS; i < splits.size() && !foundSplit; ++i )
2733                 if (( foundSplit = splits[ i ]) &&
2734                     ( n2->IsLinked( foundSplit.FirstNode()->_intPoint )))
2735                 {
2736                   iS = i - 1;
2737                 }
2738                 else
2739                 {
2740                   foundSplit._link = 0;
2741                 }
2742               if ( foundSplit )
2743               {
2744                 if ( n2 != foundSplit.FirstNode() )
2745                 {
2746                   polygon->AddPolyLink( n2, foundSplit.FirstNode() );
2747                   n2 = foundSplit.FirstNode();
2748                 }
2749                 continue;
2750               }
2751               else
2752               {
2753                 if ( n2->IsLinked( nFirst->_intPoint ))
2754                   break;
2755                 polygon->AddPolyLink( n2, n1, prevPolyg );
2756               }
2757             }
2758           } // if ( n1 != n2 )
2759
2760           polygon->_links.push_back( split );
2761           split._link = 0;
2762           --nbSplits;
2763           n2 = polygon->_links.back().LastNode();
2764
2765         } // loop on splits
2766
2767         if ( nFirst != n2 ) // close a polygon
2768         {
2769           if ( !findChain( n2, nFirst, quad, chainNodes ))
2770           {
2771             if ( !closePolygon( polygon, chainNodes ))
2772               if ( !isImplementEdges() )
2773                 chainNodes.push_back( nFirst );
2774           }
2775           for ( size_t i = 1; i < chainNodes.size(); ++i )
2776           {
2777             polygon->AddPolyLink( chainNodes[i-1], chainNodes[i], prevPolyg );
2778             nbUsedEdgeNodes += bool( chainNodes[i]->IsUsedInFace( polygon ));
2779           }
2780         }
2781
2782         if ( polygon->_links.size() < 3 && nbSplits > 0 )
2783         {
2784           polygon->_polyLinks.clear();
2785           polygon->_links.clear();
2786         }
2787       } // while ( nbSplits > 0 )
2788
2789       if ( polygon->_links.size() < 3 )
2790       {
2791         _polygons.pop_back();
2792       }
2793     }  // loop on 6 hexahedron sides
2794
2795     // Create polygons closing holes in a polyhedron
2796     // ----------------------------------------------
2797
2798     // clear _usedInFace
2799     for ( size_t iN = 0; iN < _intNodes.size(); ++iN )
2800       _intNodes[ iN ]._usedInFace = 0;
2801
2802     // add polygons to their links and mark used nodes
2803     for ( size_t iP = 0; iP < _polygons.size(); ++iP )
2804     {
2805       _Face& polygon = _polygons[ iP ];
2806       for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
2807       {
2808         polygon._links[ iL ].AddFace( &polygon );
2809         polygon._links[ iL ].FirstNode()->_usedInFace = &polygon;
2810       }
2811     }
2812     // find free links
2813     vector< _OrientedLink* > freeLinks;
2814     freeLinks.reserve(20);
2815     for ( size_t iP = 0; iP < _polygons.size(); ++iP )
2816     {
2817       _Face& polygon = _polygons[ iP ];
2818       for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
2819         if ( polygon._links[ iL ].NbFaces() < 2 )
2820           freeLinks.push_back( & polygon._links[ iL ]);
2821     }
2822     int nbFreeLinks = freeLinks.size();
2823     if ( nbFreeLinks == 1 ) return false;
2824
2825     // put not used intersection nodes to _vIntNodes
2826     int nbVertexNodes = 0; // nb not used vertex nodes
2827     {
2828       for ( size_t iN = 0; iN < _vIntNodes.size(); ++iN )
2829         nbVertexNodes += ( !_vIntNodes[ iN ]->IsUsedInFace() );
2830
2831       const double tol = 1e-3 * Min( Min( _sideLength[0], _sideLength[1] ), _sideLength[0] );
2832       for ( size_t iN = _nbFaceIntNodes; iN < _intNodes.size(); ++iN )
2833       {
2834         if ( _intNodes[ iN ].IsUsedInFace() ) continue;
2835         if ( dynamic_cast< const F_IntersectPoint* >( _intNodes[ iN ]._intPoint )) continue;
2836         _Node* equalNode =
2837           findEqualNode( _vIntNodes, _intNodes[ iN ].EdgeIntPnt(), tol*tol );
2838         if ( !equalNode )
2839         {
2840           _vIntNodes.push_back( &_intNodes[ iN ]);
2841           ++nbVertexNodes;
2842         }
2843       }
2844     }
2845
2846     set<TGeomID> usedFaceIDs;
2847     vector< TGeomID > faces;
2848     TGeomID curFace = 0;
2849     const size_t nbQuadPolygons = _polygons.size();
2850     E_IntersectPoint ipTmp;
2851
2852     // create polygons by making closed chains of free links
2853     size_t iPolygon = _polygons.size();
2854     while ( nbFreeLinks > 0 )
2855     {
2856       if ( iPolygon == _polygons.size() )
2857       {
2858         _polygons.resize( _polygons.size() + 1 );
2859         _polygons[ iPolygon ]._polyLinks.reserve( 20 );
2860         _polygons[ iPolygon ]._links.reserve( 20 );
2861       }
2862       _Face& polygon = _polygons[ iPolygon ];
2863
2864       _OrientedLink* curLink = 0;
2865       _Node*         curNode;
2866       if (( !hasEdgeIntersections ) ||
2867           ( nbFreeLinks < 4 && nbVertexNodes == 0 ))
2868       {
2869         // get a remaining link to start from
2870         for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2871           if (( curLink = freeLinks[ iL ] ))
2872             freeLinks[ iL ] = 0;
2873         polygon._links.push_back( *curLink );
2874         --nbFreeLinks;
2875         do
2876         {
2877           // find all links connected to curLink
2878           curNode = curLink->FirstNode();
2879           curLink = 0;
2880           for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2881             if ( freeLinks[ iL ] && freeLinks[ iL ]->LastNode() == curNode )
2882             {
2883               curLink = freeLinks[ iL ];
2884               freeLinks[ iL ] = 0;
2885               --nbFreeLinks;
2886               polygon._links.push_back( *curLink );
2887             }
2888         } while ( curLink );
2889       }
2890       else // there are intersections with EDGEs
2891       {
2892         // get a remaining link to start from, one lying on minimal nb of FACEs
2893         {
2894           typedef pair< TGeomID, int > TFaceOfLink;
2895           TFaceOfLink faceOfLink( -1, -1 );
2896           TFaceOfLink facesOfLink[3] = { faceOfLink, faceOfLink, faceOfLink };
2897           for ( size_t iL = 0; iL < freeLinks.size(); ++iL )
2898             if ( freeLinks[ iL ] )
2899             {
2900               faces = freeLinks[ iL ]->GetNotUsedFace( usedFaceIDs );
2901               if ( faces.size() == 1 )
2902               {
2903                 faceOfLink = TFaceOfLink( faces[0], iL );
2904                 if ( !freeLinks[ iL ]->HasEdgeNodes() )
2905                   break;
2906                 facesOfLink[0] = faceOfLink;
2907               }
2908               else if ( facesOfLink[0].first < 0 )
2909               {
2910                 faceOfLink = TFaceOfLink(( faces.empty() ? -1 : faces[0]), iL );
2911                 facesOfLink[ 1 + faces.empty() ] = faceOfLink;
2912               }
2913             }
2914           for ( int i = 0; faceOfLink.first < 0 && i < 3; ++i )
2915             faceOfLink = facesOfLink[i];
2916
2917           if ( faceOfLink.first < 0 ) // all faces used
2918           {
2919             for ( size_t iL = 0; iL < freeLinks.size() && faceOfLink.first < 1; ++iL )
2920               if (( curLink = freeLinks[ iL ]))
2921               {
2922                 faceOfLink.first = 
2923                   curLink->FirstNode()->IsLinked( curLink->LastNode()->_intPoint );
2924                 faceOfLink.second = iL;
2925               }
2926             usedFaceIDs.clear();
2927           }
2928           curFace = faceOfLink.first;
2929           curLink = freeLinks[ faceOfLink.second ];
2930           freeLinks[ faceOfLink.second ] = 0;
2931         }
2932         usedFaceIDs.insert( curFace );
2933         polygon._links.push_back( *curLink );
2934         --nbFreeLinks;
2935
2936         // find all links lying on a curFace
2937         do
2938         {
2939           // go forward from curLink
2940           curNode = curLink->LastNode();
2941           curLink = 0;
2942           for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2943             if ( freeLinks[ iL ] &&
2944                  freeLinks[ iL ]->FirstNode() == curNode &&
2945                  freeLinks[ iL ]->LastNode()->IsOnFace( curFace ))
2946             {
2947               curLink = freeLinks[ iL ];
2948               freeLinks[ iL ] = 0;
2949               polygon._links.push_back( *curLink );
2950               --nbFreeLinks;
2951             }
2952         } while ( curLink );
2953
2954         std::reverse( polygon._links.begin(), polygon._links.end() );
2955
2956         curLink = & polygon._links.back();
2957         do
2958         {
2959           // go backward from curLink
2960           curNode = curLink->FirstNode();
2961           curLink = 0;
2962           for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2963             if ( freeLinks[ iL ] &&
2964                  freeLinks[ iL ]->LastNode() == curNode &&
2965                  freeLinks[ iL ]->FirstNode()->IsOnFace( curFace ))
2966             {
2967               curLink = freeLinks[ iL ];
2968               freeLinks[ iL ] = 0;
2969               polygon._links.push_back( *curLink );
2970               --nbFreeLinks;
2971             }
2972         } while ( curLink );
2973
2974         curNode = polygon._links.back().FirstNode();
2975
2976         if ( polygon._links[0].LastNode() != curNode )
2977         {
2978           if ( nbVertexNodes > 0 )
2979           {
2980             // add links with _vIntNodes if not already used
2981             chainNodes.clear();
2982             for ( size_t iN = 0; iN < _vIntNodes.size(); ++iN )
2983               if ( !_vIntNodes[ iN ]->IsUsedInFace() &&
2984                    _vIntNodes[ iN ]->IsOnFace( curFace ))
2985               {
2986                 _vIntNodes[ iN ]->_usedInFace = &polygon;
2987                 chainNodes.push_back( _vIntNodes[ iN ] );
2988               }
2989             if ( chainNodes.size() > 1 &&
2990                  curFace != _grid->PseudoIntExtFaceID() ) /////// TODO
2991             {
2992               sortVertexNodes( chainNodes, curNode, curFace );
2993             }
2994             for ( size_t i = 0; i < chainNodes.size(); ++i )
2995             {
2996               polygon.AddPolyLink( chainNodes[ i ], curNode );
2997               curNode = chainNodes[ i ];
2998               freeLinks.push_back( &polygon._links.back() );
2999               ++nbFreeLinks;
3000             }
3001             nbVertexNodes -= chainNodes.size();
3002           }
3003           // if ( polygon._links.size() > 1 )
3004           {
3005             polygon.AddPolyLink( polygon._links[0].LastNode(), curNode );
3006             freeLinks.push_back( &polygon._links.back() );
3007             ++nbFreeLinks;
3008           }
3009         }
3010       } // if there are intersections with EDGEs
3011
3012       if ( polygon._links.size() < 2 ||
3013            polygon._links[0].LastNode() != polygon._links.back().FirstNode() )
3014         return false; // closed polygon not found -> invalid polyhedron
3015
3016       if ( polygon._links.size() == 2 )
3017       {
3018         if ( freeLinks.back() == &polygon._links.back() )
3019         {
3020           freeLinks.pop_back();
3021           --nbFreeLinks;
3022         }
3023         if ( polygon._links.front().NbFaces() > 0 )
3024           polygon._links.back().AddFace( polygon._links.front()._link->_faces[0] );
3025         if ( polygon._links.back().NbFaces() > 0 )
3026           polygon._links.front().AddFace( polygon._links.back()._link->_faces[0] );
3027
3028         if ( iPolygon == _polygons.size()-1 )
3029           _polygons.pop_back();
3030       }
3031       else // polygon._links.size() >= 2
3032       {
3033         // add polygon to its links
3034         for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
3035         {
3036           polygon._links[ iL ].AddFace( &polygon );
3037           polygon._links[ iL ].Reverse();
3038         }
3039         if ( /*hasEdgeIntersections &&*/ iPolygon == _polygons.size() - 1 )
3040         {
3041           // check that a polygon does not lie on a hexa side
3042           coplanarPolyg = 0;
3043           for ( size_t iL = 0; iL < polygon._links.size() && !coplanarPolyg; ++iL )
3044           {
3045             if ( polygon._links[ iL ].NbFaces() < 2 )
3046               continue; // it's a just added free link
3047             // look for a polygon made on a hexa side and sharing
3048             // two or more haxa links
3049             size_t iL2;
3050             coplanarPolyg = polygon._links[ iL ]._link->_faces[0];
3051             for ( iL2 = iL + 1; iL2 < polygon._links.size(); ++iL2 )
3052               if ( polygon._links[ iL2 ]._link->_faces[0] == coplanarPolyg &&
3053                    !coplanarPolyg->IsPolyLink( polygon._links[ iL  ]) &&
3054                    !coplanarPolyg->IsPolyLink( polygon._links[ iL2 ]) &&
3055                    coplanarPolyg < & _polygons[ nbQuadPolygons ])
3056                 break;
3057             if ( iL2 == polygon._links.size() )
3058               coplanarPolyg = 0;
3059           }
3060           if ( coplanarPolyg ) // coplanar polygon found
3061           {
3062             freeLinks.resize( freeLinks.size() - polygon._polyLinks.size() );
3063             nbFreeLinks -= polygon._polyLinks.size();
3064
3065             // an E_IntersectPoint used to mark nodes of coplanarPolyg
3066             // as lying on curFace while they are not at intersection with geometry
3067             ipTmp._faceIDs.resize(1);
3068             ipTmp._faceIDs[0] = curFace;
3069
3070             // fill freeLinks with links not shared by coplanarPolyg and polygon
3071             for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
3072               if ( polygon._links[ iL ]._link->_faces[1] &&
3073                    polygon._links[ iL ]._link->_faces[0] != coplanarPolyg )
3074               {
3075                 _Face* p = polygon._links[ iL ]._link->_faces[0];
3076                 for ( size_t iL2 = 0; iL2 < p->_links.size(); ++iL2 )
3077                   if ( p->_links[ iL2 ]._link == polygon._links[ iL ]._link )
3078                   {
3079                     freeLinks.push_back( & p->_links[ iL2 ] );
3080                     ++nbFreeLinks;
3081                     freeLinks.back()->RemoveFace( &polygon );
3082                     break;
3083                   }
3084               }
3085             for ( size_t iL = 0; iL < coplanarPolyg->_links.size(); ++iL )
3086               if ( coplanarPolyg->_links[ iL ]._link->_faces[1] &&
3087                    coplanarPolyg->_links[ iL ]._link->_faces[1] != &polygon )
3088               {
3089                 _Face* p = coplanarPolyg->_links[ iL ]._link->_faces[0];
3090                 if ( p == coplanarPolyg )
3091                   p = coplanarPolyg->_links[ iL ]._link->_faces[1];
3092                 for ( size_t iL2 = 0; iL2 < p->_links.size(); ++iL2 )
3093                   if ( p->_links[ iL2 ]._link == coplanarPolyg->_links[ iL ]._link )
3094                   {
3095                     // set links of coplanarPolyg in place of used freeLinks
3096                     // to re-create coplanarPolyg next
3097                     size_t iL3 = 0;
3098                     for ( ; iL3 < freeLinks.size() && freeLinks[ iL3 ]; ++iL3 );
3099                     if ( iL3 < freeLinks.size() )
3100                       freeLinks[ iL3 ] = ( & p->_links[ iL2 ] );
3101                     else
3102                       freeLinks.push_back( & p->_links[ iL2 ] );
3103                     ++nbFreeLinks;
3104                     freeLinks[ iL3 ]->RemoveFace( coplanarPolyg );
3105                     //  mark nodes of coplanarPolyg as lying on curFace
3106                     for ( int iN = 0; iN < 2; ++iN )
3107                     {
3108                       _Node* n = freeLinks[ iL3 ]->_link->_nodes[ iN ];
3109                       if ( n->_intPoint ) n->_intPoint->Add( ipTmp._faceIDs );
3110                       else                n->_intPoint = &ipTmp;
3111                     }
3112                     break;
3113                   }
3114               }
3115             // set coplanarPolyg to be re-created next
3116             for ( size_t iP = 0; iP < _polygons.size(); ++iP )
3117               if ( coplanarPolyg == & _polygons[ iP ] )
3118               {
3119                 iPolygon = iP;
3120                 _polygons[ iPolygon ]._links.clear();
3121                 _polygons[ iPolygon ]._polyLinks.clear();
3122                 break;
3123               }
3124             _polygons.pop_back();
3125             usedFaceIDs.erase( curFace );
3126             continue;
3127           } // if ( coplanarPolyg )
3128         } // if ( hasEdgeIntersections ) - search for coplanarPolyg
3129
3130         iPolygon = _polygons.size();
3131
3132       } // end of case ( polygon._links.size() > 2 )
3133     } // while ( nbFreeLinks > 0 )
3134
3135     // check volume size
3136     _hasTooSmall = ! checkPolyhedronSize( intFlag & IS_CUT_BY_INTERNAL_FACE );
3137
3138     for ( size_t i = 0; i < 8; ++i )
3139       if ( _hexNodes[ i ]._intPoint == &ipTmp )
3140         _hexNodes[ i ]._intPoint = 0;
3141
3142     if ( _hasTooSmall )
3143       return false; // too small volume
3144
3145     // create a classic cell if possible
3146
3147     int nbPolygons = 0;
3148     for ( size_t iF = 0; iF < _polygons.size(); ++iF )
3149       nbPolygons += (_polygons[ iF ]._links.size() > 0 );
3150
3151     //const int nbNodes = _nbCornerNodes + nbIntersections;
3152     int nbNodes = 0;
3153     for ( size_t i = 0; i < 8; ++i )
3154       nbNodes += _hexNodes[ i ].IsUsedInFace();
3155     for ( size_t i = 0; i < _intNodes.size(); ++i )
3156       nbNodes += _intNodes[ i ].IsUsedInFace();
3157
3158     bool isClassicElem = false;
3159     if (      nbNodes == 8 && nbPolygons == 6 ) isClassicElem = addHexa();
3160     else if ( nbNodes == 4 && nbPolygons == 4 ) isClassicElem = addTetra();
3161     else if ( nbNodes == 6 && nbPolygons == 5 ) isClassicElem = addPenta();
3162     else if ( nbNodes == 5 && nbPolygons == 5 ) isClassicElem = addPyra ();
3163     if ( !isClassicElem )
3164     {
3165       _volumeDefs._nodes.clear();
3166       _volumeDefs._quantities.clear();
3167
3168       for ( size_t iF = 0; iF < _polygons.size(); ++iF )
3169       {
3170         const size_t nbLinks = _polygons[ iF ]._links.size();
3171         if ( nbLinks == 0 ) continue;
3172         _volumeDefs._quantities.push_back( nbLinks );
3173         for ( size_t iL = 0; iL < nbLinks; ++iL )
3174           _volumeDefs._nodes.push_back( _polygons[ iF ]._links[ iL ].FirstNode() );
3175       }
3176     }
3177     _volumeDefs._solidID = solid->ID();
3178
3179     return !_volumeDefs._nodes.empty();
3180   }
3181   //================================================================================
3182   /*!
3183    * \brief Create elements in the mesh
3184    */
3185   int Hexahedron::MakeElements(SMESH_MesherHelper&                      helper,
3186                                const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap)
3187   {
3188     SMESHDS_Mesh* mesh = helper.GetMeshDS();
3189
3190     CellsAroundLink c( _grid, 0 );
3191     const size_t nbGridCells = c._nbCells[0] * c._nbCells[1] * c._nbCells[2];
3192     vector< Hexahedron* > allHexa( nbGridCells, 0 );
3193     int nbIntHex = 0;
3194
3195     // set intersection nodes from GridLine's to links of allHexa
3196     int i,j,k, cellIndex;
3197     for ( int iDir = 0; iDir < 3; ++iDir )
3198     {
3199       // loop on GridLine's parallel to iDir
3200       LineIndexer lineInd = _grid->GetLineIndexer( iDir );
3201       CellsAroundLink fourCells( _grid, iDir );
3202       for ( ; lineInd.More(); ++lineInd )
3203       {
3204         GridLine& line = _grid->_lines[ iDir ][ lineInd.LineIndex() ];
3205         multiset< F_IntersectPoint >::const_iterator ip = line._intPoints.begin();
3206         for ( ; ip != line._intPoints.end(); ++ip )
3207         {
3208           // if ( !ip->_node ) continue; // intersection at a grid node
3209           lineInd.SetIndexOnLine( ip->_indexOnLine );
3210           fourCells.Init( lineInd.I(), lineInd.J(), lineInd.K() );
3211           for ( int iL = 0; iL < 4; ++iL ) // loop on 4 cells sharing a link
3212           {
3213             if ( !fourCells.GetCell( iL, i,j,k, cellIndex ))
3214               continue;
3215             Hexahedron *& hex = allHexa[ cellIndex ];
3216             if ( !hex)
3217             {
3218               hex = new Hexahedron( *this, i, j, k, cellIndex );
3219               ++nbIntHex;
3220             }
3221             const int iLink = iL + iDir * 4;
3222             hex->_hexLinks[iLink]._fIntPoints.push_back( &(*ip) );
3223             hex->_nbFaceIntNodes += bool( ip->_node );
3224           }
3225         }
3226       }
3227     }
3228
3229     // implement geom edges into the mesh
3230     addEdges( helper, allHexa, edge2faceIDsMap );
3231
3232     // add not split hexahedra to the mesh
3233     int nbAdded = 0;
3234     vector< Hexahedron* > intHexa; intHexa.reserve( nbIntHex );
3235     vector< const SMDS_MeshElement* > boundaryVolumes; boundaryVolumes.reserve( nbIntHex * 1.1 );
3236     for ( size_t i = 0; i < allHexa.size(); ++i )
3237     {
3238       // initialize this by not cut allHexa[ i ]
3239       Hexahedron * & hex = allHexa[ i ];
3240       if ( hex ) // split hexahedron
3241       {
3242         intHexa.push_back( hex );
3243         if ( hex->_nbFaceIntNodes > 0 || hex->_eIntPoints.size() > 0 )
3244           continue; // treat intersected hex later in parallel
3245         this->init( hex->_i, hex->_j, hex->_k );
3246       }
3247       else
3248       {
3249         this->init( i ); // == init(i,j,k)
3250       }
3251       if (( _nbCornerNodes == 8 ) &&
3252           ( _nbBndNodes < _nbCornerNodes || !isInHole() ))
3253       {
3254         // order of _hexNodes is defined by enum SMESH_Block::TShapeID
3255         SMDS_MeshElement* el =
3256           mesh->AddVolume( _hexNodes[0].Node(), _hexNodes[2].Node(),
3257                            _hexNodes[3].Node(), _hexNodes[1].Node(),
3258                            _hexNodes[4].Node(), _hexNodes[6].Node(),
3259                            _hexNodes[7].Node(), _hexNodes[5].Node() );
3260         TGeomID solidID = 0;
3261         if ( _nbBndNodes < _nbCornerNodes )
3262         {
3263           for ( int iN = 0; iN < 8 &&  !solidID; ++iN )
3264             if ( !_hexNodes[iN]._intPoint ) // no intersection
3265               solidID = _hexNodes[iN].Node()->GetShapeID();
3266         }
3267         else
3268         {
3269           solidID = getSolids()[0];
3270         }
3271         mesh->SetMeshElementOnShape( el, solidID );
3272         ++nbAdded;
3273         if ( hex )
3274           intHexa.pop_back();
3275         if ( _grid->_toCreateFaces && _nbBndNodes >= 3 )
3276         {
3277           boundaryVolumes.push_back( el );
3278           el->setIsMarked( true );
3279         }
3280       }
3281       else if ( _nbCornerNodes > 3 && !hex )
3282       {
3283         // all intersection of hex with geometry are at grid nodes
3284         hex = new Hexahedron( *this, _i, _j, _k, i );
3285         intHexa.push_back( hex );
3286       }
3287     }
3288
3289     // add elements resulted from hexadron intersection
3290 #ifdef WITH_TBB
3291     tbb::parallel_for ( tbb::blocked_range<size_t>( 0, intHexa.size() ),
3292                         ParallelHexahedron( intHexa ),
3293                         tbb::simple_partitioner()); // ComputeElements() is called here
3294     for ( size_t i = 0; i < intHexa.size(); ++i )
3295       if ( Hexahedron * hex = intHexa[ i ] )
3296         nbAdded += hex->addVolumes( helper );
3297 #else
3298     for ( size_t i = 0; i < intHexa.size(); ++i )
3299       if ( Hexahedron * hex = intHexa[ i ] )
3300       {
3301         hex->ComputeElements();
3302         nbAdded += hex->addVolumes( helper );
3303       }
3304 #endif
3305
3306     // fill boundaryVolumes with volumes neighboring too small skipped volumes
3307     if ( _grid->_toCreateFaces )
3308     {
3309       for ( size_t i = 0; i < intHexa.size(); ++i )
3310         if ( Hexahedron * hex = intHexa[ i ] )
3311           hex->getBoundaryElems( boundaryVolumes );
3312     }
3313
3314     // create boundary mesh faces
3315     addFaces( helper, boundaryVolumes );
3316
3317     // create mesh edges
3318     addSegments( helper, edge2faceIDsMap );
3319
3320     for ( size_t i = 0; i < allHexa.size(); ++i )
3321       if ( allHexa[ i ] )
3322         delete allHexa[ i ];
3323
3324     return nbAdded;
3325   }
3326
3327   //================================================================================
3328   /*!
3329    * \brief Implements geom edges into the mesh
3330    */
3331   void Hexahedron::addEdges(SMESH_MesherHelper&                      helper,
3332                             vector< Hexahedron* >&                   hexes,
3333                             const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap)
3334   {
3335     if ( edge2faceIDsMap.empty() ) return;
3336
3337     // Prepare planes for intersecting with EDGEs
3338     GridPlanes pln[3];
3339     {
3340       for ( int iDirZ = 0; iDirZ < 3; ++iDirZ ) // iDirZ gives normal direction to planes
3341       {
3342         GridPlanes& planes = pln[ iDirZ ];
3343         int iDirX = ( iDirZ + 1 ) % 3;
3344         int iDirY = ( iDirZ + 2 ) % 3;
3345         planes._zNorm  = ( _grid->_axes[ iDirX ] ^ _grid->_axes[ iDirY ] ).Normalized();
3346         planes._zProjs.resize ( _grid->_coords[ iDirZ ].size() );
3347         planes._zProjs [0] = 0;
3348         const double       zFactor = _grid->_axes[ iDirZ ] * planes._zNorm;
3349         const vector< double > & u = _grid->_coords[ iDirZ ];
3350         for ( size_t i = 1; i < planes._zProjs.size(); ++i )
3351         {
3352           planes._zProjs [i] = zFactor * ( u[i] - u[0] );
3353         }
3354       }
3355     }
3356     const double deflection = _grid->_minCellSize / 20.;
3357     const double tol        = _grid->_tol;
3358     E_IntersectPoint ip;
3359
3360     TColStd_MapOfInteger intEdgeIDs; // IDs of not shared INTERNAL EDGES
3361
3362     // Intersect EDGEs with the planes
3363     map< TGeomID, vector< TGeomID > >::const_iterator e2fIt = edge2faceIDsMap.begin();
3364     for ( ; e2fIt != edge2faceIDsMap.end(); ++e2fIt )
3365     {
3366       const TGeomID  edgeID = e2fIt->first;
3367       const TopoDS_Edge & E = TopoDS::Edge( _grid->Shape( edgeID ));
3368       BRepAdaptor_Curve curve( E );
3369       TopoDS_Vertex v1 = helper.IthVertex( 0, E, false );
3370       TopoDS_Vertex v2 = helper.IthVertex( 1, E, false );
3371
3372       ip._faceIDs = e2fIt->second;
3373       ip._shapeID = edgeID;
3374
3375       bool isInternal = ( ip._faceIDs.size() == 1 && _grid->IsInternal( edgeID ));
3376       if ( isInternal )
3377       {
3378         intEdgeIDs.Add( edgeID );
3379         intEdgeIDs.Add( _grid->ShapeID( v1 ));
3380         intEdgeIDs.Add( _grid->ShapeID( v2 ));
3381       }
3382
3383       // discretize the EDGE
3384       GCPnts_UniformDeflection discret( curve, deflection, true );
3385       if ( !discret.IsDone() || discret.NbPoints() < 2 )
3386         continue;
3387
3388       // perform intersection
3389       E_IntersectPoint* eip, *vip;
3390       for ( int iDirZ = 0; iDirZ < 3; ++iDirZ )
3391       {
3392         GridPlanes& planes = pln[ iDirZ ];
3393         int      iDirX = ( iDirZ + 1 ) % 3;
3394         int      iDirY = ( iDirZ + 2 ) % 3;
3395         double    xLen = _grid->_coords[ iDirX ].back() - _grid->_coords[ iDirX ][0];
3396         double    yLen = _grid->_coords[ iDirY ].back() - _grid->_coords[ iDirY ][0];
3397         double    zLen = _grid->_coords[ iDirZ ].back() - _grid->_coords[ iDirZ ][0];
3398         int dIJK[3], d000[3] = { 0,0,0 };
3399         double o[3] = { _grid->_coords[0][0],
3400                         _grid->_coords[1][0],
3401                         _grid->_coords[2][0] };
3402
3403         // locate the 1st point of a segment within the grid
3404         gp_XYZ p1     = discret.Value( 1 ).XYZ();
3405         double u1     = discret.Parameter( 1 );
3406         double zProj1 = planes._zNorm * ( p1 - _grid->_origin );
3407
3408         _grid->ComputeUVW( p1, ip._uvw );
3409         int iX1 = int(( ip._uvw[iDirX] - o[iDirX]) / xLen * (_grid->_coords[ iDirX ].size() - 1));
3410         int iY1 = int(( ip._uvw[iDirY] - o[iDirY]) / yLen * (_grid->_coords[ iDirY ].size() - 1));
3411         int iZ1 = int(( ip._uvw[iDirZ] - o[iDirZ]) / zLen * (_grid->_coords[ iDirZ ].size() - 1));
3412         locateValue( iX1, ip._uvw[iDirX], _grid->_coords[ iDirX ], dIJK[ iDirX ], tol );
3413         locateValue( iY1, ip._uvw[iDirY], _grid->_coords[ iDirY ], dIJK[ iDirY ], tol );
3414         locateValue( iZ1, ip._uvw[iDirZ], _grid->_coords[ iDirZ ], dIJK[ iDirZ ], tol );
3415
3416         int ijk[3]; // grid index where a segment intersects a plane
3417         ijk[ iDirX ] = iX1;
3418         ijk[ iDirY ] = iY1;
3419         ijk[ iDirZ ] = iZ1;
3420
3421         // add the 1st vertex point to a hexahedron
3422         if ( iDirZ == 0 )
3423         {
3424           ip._point   = p1;
3425           ip._shapeID = _grid->ShapeID( v1 );
3426           vip = _grid->Add( ip );
3427           if ( isInternal )
3428             vip->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
3429           if ( !addIntersection( vip, hexes, ijk, d000 ))
3430             _grid->Remove( vip );
3431           ip._shapeID = edgeID;
3432         }
3433         for ( int iP = 2; iP <= discret.NbPoints(); ++iP )
3434         {
3435           // locate the 2nd point of a segment within the grid
3436           gp_XYZ p2     = discret.Value( iP ).XYZ();
3437           double u2     = discret.Parameter( iP );
3438           double zProj2 = planes._zNorm * ( p2 - _grid->_origin );
3439           int    iZ2    = iZ1;
3440           if ( Abs( zProj2 - zProj1 ) > std::numeric_limits<double>::min() )
3441           {
3442             locateValue( iZ2, zProj2, planes._zProjs, dIJK[ iDirZ ], tol );
3443
3444             // treat intersections with planes between 2 end points of a segment
3445             int dZ = ( iZ1 <= iZ2 ) ? +1 : -1;
3446             int iZ = iZ1 + ( iZ1 < iZ2 );
3447             for ( int i = 0, nb = Abs( iZ1 - iZ2 ); i < nb; ++i, iZ += dZ )
3448             {
3449               ip._point = findIntPoint( u1, zProj1, u2, zProj2,
3450                                         planes._zProjs[ iZ ],
3451                                         curve, planes._zNorm, _grid->_origin );
3452               _grid->ComputeUVW( ip._point.XYZ(), ip._uvw );
3453               locateValue( ijk[iDirX], ip._uvw[iDirX], _grid->_coords[iDirX], dIJK[iDirX], tol );
3454               locateValue( ijk[iDirY], ip._uvw[iDirY], _grid->_coords[iDirY], dIJK[iDirY], tol );
3455               ijk[ iDirZ ] = iZ;
3456
3457               // add ip to hex "above" the plane
3458               eip = _grid->Add( ip );
3459               if ( isInternal )
3460                 eip->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
3461               dIJK[ iDirZ ] = 0;
3462               bool added = addIntersection( eip, hexes, ijk, dIJK);
3463
3464               // add ip to hex "below" the plane
3465               ijk[ iDirZ ] = iZ-1;
3466               if ( !addIntersection( eip, hexes, ijk, dIJK ) &&
3467                    !added )
3468                 _grid->Remove( eip );
3469             }
3470           }
3471           iZ1    = iZ2;
3472           p1     = p2;
3473           u1     = u2;
3474           zProj1 = zProj2;
3475         }
3476         // add the 2nd vertex point to a hexahedron
3477         if ( iDirZ == 0 )
3478         {
3479           ip._point   = p1;
3480           ip._shapeID = _grid->ShapeID( v2 );
3481           _grid->ComputeUVW( p1, ip._uvw );
3482           locateValue( ijk[iDirX], ip._uvw[iDirX], _grid->_coords[iDirX], dIJK[iDirX], tol );
3483           locateValue( ijk[iDirY], ip._uvw[iDirY], _grid->_coords[iDirY], dIJK[iDirY], tol );
3484           ijk[ iDirZ ] = iZ1;
3485           bool sameV = ( v1.IsSame( v2 ));
3486           if ( !sameV )
3487             vip = _grid->Add( ip );
3488           if ( isInternal && !sameV )
3489             vip->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
3490           if ( !addIntersection( vip, hexes, ijk, d000 ) && !sameV )
3491             _grid->Remove( vip );
3492           ip._shapeID = edgeID;
3493         }
3494       } // loop on 3 grid directions
3495     } // loop on EDGEs
3496
3497
3498     if ( intEdgeIDs.Size() > 0 )
3499       cutByExtendedInternal( hexes, intEdgeIDs );
3500
3501     return;
3502   }
3503
3504   //================================================================================
3505   /*!
3506    * \brief Fully cut hexes that are partially cut by INTERNAL FACE.
3507    *        Cut them by extended INTERNAL FACE.
3508    */
3509   void Hexahedron::cutByExtendedInternal( std::vector< Hexahedron* >& hexes,
3510                                           const TColStd_MapOfInteger& intEdgeIDs )
3511   {
3512     IntAna_IntConicQuad intersection;
3513     SMESHDS_Mesh* meshDS = _grid->_helper->GetMeshDS();
3514     const double tol2 = _grid->_tol * _grid->_tol;
3515
3516     for ( size_t iH = 0; iH < hexes.size(); ++iH )
3517     {
3518       Hexahedron* hex = hexes[ iH ];
3519       if ( !hex || hex->_eIntPoints.size() < 2 )
3520         continue;
3521       if ( !intEdgeIDs.Contains( hex->_eIntPoints.back()->_shapeID ))
3522         continue;
3523
3524       // get 3 points on INTERNAL FACE to construct a cutting plane
3525       gp_Pnt p1 = hex->_eIntPoints[0]->_point;
3526       gp_Pnt p2 = hex->_eIntPoints[1]->_point;
3527       gp_Pnt p3 = hex->mostDistantInternalPnt( iH, p1, p2 );
3528
3529       gp_Vec norm = gp_Vec( p1, p2 ) ^ gp_Vec( p1, p3 );
3530       gp_Pln pln;
3531       try {
3532         pln = gp_Pln( p1, norm );
3533       }
3534       catch(...)
3535       {
3536         continue;
3537       }
3538
3539       TGeomID intFaceID = hex->_eIntPoints.back()->_faceIDs.front(); // FACE being "extended"
3540       TGeomID   solidID = _grid->GetSolid( intFaceID )->ID();
3541
3542       // cut links by the plane
3543       //bool isCut = false;
3544       for ( int iLink = 0; iLink < 12; ++iLink )
3545       {
3546         _Link& link = hex->_hexLinks[ iLink ];
3547         if ( !link._fIntPoints.empty() )
3548         {
3549           // if ( link._fIntPoints[0]->_faceIDs.back() == _grid->PseudoIntExtFaceID() )
3550           //   isCut = true;
3551           continue; // already cut link
3552         }
3553         if ( !link._nodes[0]->Node() ||
3554              !link._nodes[1]->Node() )
3555           continue; // outside link
3556
3557         if ( link._nodes[0]->IsOnFace( intFaceID ))
3558         {
3559           if ( link._nodes[0]->_intPoint->_faceIDs.back() != _grid->PseudoIntExtFaceID() )
3560             if ( p1.SquareDistance( link._nodes[0]->Point() ) < tol2  ||
3561                  p2.SquareDistance( link._nodes[0]->Point() ) < tol2 )
3562               link._nodes[0]->_intPoint->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
3563           continue; // link is cut by FACE being "extended"
3564         }
3565         if ( link._nodes[1]->IsOnFace( intFaceID ))
3566         {
3567           if ( link._nodes[1]->_intPoint->_faceIDs.back() != _grid->PseudoIntExtFaceID() )
3568             if ( p1.SquareDistance( link._nodes[1]->Point() ) < tol2  ||
3569                  p2.SquareDistance( link._nodes[1]->Point() ) < tol2 )
3570               link._nodes[1]->_intPoint->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
3571           continue; // link is cut by FACE being "extended"
3572         }
3573         gp_Pnt p4 = link._nodes[0]->Point();
3574         gp_Pnt p5 = link._nodes[1]->Point();
3575         gp_Lin line( p4, gp_Vec( p4, p5 ));
3576
3577         intersection.Perform( line, pln );
3578         if ( !intersection.IsDone() ||
3579              intersection.IsInQuadric() ||
3580              intersection.IsParallel() ||
3581              intersection.NbPoints() < 1 )
3582           continue;
3583
3584         double u = intersection.ParamOnConic(1);
3585         if ( u + _grid->_tol < 0 )
3586           continue;
3587         int       iDir = iLink / 4;
3588         int      index = (&hex->_i)[iDir];
3589         double linkLen = _grid->_coords[iDir][index+1] - _grid->_coords[iDir][index];
3590         if ( u - _grid->_tol > linkLen )
3591           continue;
3592
3593         if ( u < _grid->_tol ||
3594              u > linkLen - _grid->_tol ) // intersection at grid node
3595         {
3596           int  i = ! ( u < _grid->_tol ); // [0,1]
3597           int iN = link._nodes[ i ] - hex->_hexNodes; // [0-7]
3598
3599           const F_IntersectPoint * & ip = _grid->_gridIntP[ hex->_origNodeInd + _nodeShift[iN] ];
3600           if ( !ip )
3601           {
3602             ip = _grid->_extIntPool.getNew();
3603             ip->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
3604             //ip->_transition = Trans_INTERNAL;
3605           }
3606           else if ( ip->_faceIDs.back() != _grid->PseudoIntExtFaceID() )
3607           {
3608             ip->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
3609           }
3610           hex->_nbFaceIntNodes++;
3611           //isCut = true;
3612         }
3613         else
3614         {
3615           const gp_Pnt&      p = intersection.Point( 1 );
3616           F_IntersectPoint* ip = _grid->_extIntPool.getNew();
3617           ip->_node = meshDS->AddNode( p.X(), p.Y(), p.Z() );
3618           ip->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
3619           ip->_transition = Trans_INTERNAL;
3620           meshDS->SetNodeInVolume( ip->_node, solidID );
3621
3622           CellsAroundLink fourCells( _grid, iDir );
3623           fourCells.Init( hex->_i, hex->_j, hex->_k, iLink );
3624           int i,j,k, cellIndex;
3625           for ( int iC = 0; iC < 4; ++iC ) // loop on 4 cells sharing the link
3626           {
3627             if ( !fourCells.GetCell( iC, i,j,k, cellIndex ))
3628               continue;
3629             Hexahedron * h = hexes[ cellIndex ];
3630             if ( !h )
3631               h = hexes[ cellIndex ] = new Hexahedron( *this, i, j, k, cellIndex );
3632             const int iL = iC + iDir * 4;
3633             h->_hexLinks[iL]._fIntPoints.push_back( ip );
3634             h->_nbFaceIntNodes++;
3635             //isCut = true;
3636           }
3637         }
3638       }
3639
3640       // if ( isCut )
3641       //   for ( size_t i = 0; i < hex->_eIntPoints.size(); ++i )
3642       //   {
3643       //     if ( _grid->IsInternal( hex->_eIntPoints[i]->_shapeID ) &&
3644       //          ! hex->_eIntPoints[i]->IsOnFace( _grid->PseudoIntExtFaceID() ))
3645       //       hex->_eIntPoints[i]->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
3646       //   }
3647       continue;
3648
3649     } // loop on all hexes
3650     return;
3651   }
3652
3653   //================================================================================
3654   /*!
3655    * \brief Return intersection point on INTERNAL FACE most distant from given ones
3656    */
3657   gp_Pnt Hexahedron::mostDistantInternalPnt( int hexIndex, const gp_Pnt& p1, const gp_Pnt& p2 )
3658   {
3659     gp_Pnt resultPnt = p1;
3660
3661     double maxDist2 = 0;
3662     for ( int iLink = 0; iLink < 12; ++iLink ) // check links
3663     {
3664       _Link& link = _hexLinks[ iLink ];
3665       for ( size_t i = 0; i < link._fIntPoints.size(); ++i )
3666         if ( _grid->PseudoIntExtFaceID() != link._fIntPoints[i]->_faceIDs[0] &&
3667              _grid->IsInternal( link._fIntPoints[i]->_faceIDs[0] ) &&
3668              link._fIntPoints[i]->_node )
3669         {
3670           gp_Pnt p = SMESH_NodeXYZ( link._fIntPoints[i]->_node );
3671           double d = p1.SquareDistance( p );
3672           if ( d > maxDist2 )
3673           {
3674             resultPnt = p;
3675             maxDist2  = d;
3676           }
3677           else
3678           {
3679             d = p2.SquareDistance( p );
3680             if ( d > maxDist2 )
3681             {
3682               resultPnt = p;
3683               maxDist2  = d;
3684             }
3685           }
3686         }
3687     }
3688     setIJK( hexIndex );
3689     _origNodeInd = _grid->NodeIndex( _i,_j,_k );
3690
3691     for ( size_t iN = 0; iN < 8; ++iN ) // check corners
3692     {
3693       _hexNodes[iN]._node     = _grid->_nodes   [ _origNodeInd + _nodeShift[iN] ];
3694       _hexNodes[iN]._intPoint = _grid->_gridIntP[ _origNodeInd + _nodeShift[iN] ];
3695       if ( _hexNodes[iN]._intPoint )
3696         for ( size_t iF = 0; iF < _hexNodes[iN]._intPoint->_faceIDs.size(); ++iF )
3697         {
3698           if ( _grid->IsInternal( _hexNodes[iN]._intPoint->_faceIDs[iF]))
3699           {
3700             gp_Pnt p = SMESH_NodeXYZ( _hexNodes[iN]._node );
3701             double d = p1.SquareDistance( p );
3702             if ( d > maxDist2 )
3703             {
3704               resultPnt = p;
3705               maxDist2  = d;
3706             }
3707             else
3708             {
3709               d = p2.SquareDistance( p );
3710               if ( d > maxDist2 )
3711               {
3712                 resultPnt = p;
3713                 maxDist2  = d;
3714               }
3715             }
3716           }
3717         }
3718     }
3719     if ( maxDist2 < _grid->_tol * _grid->_tol )
3720       return p1;
3721
3722     return resultPnt;
3723   }
3724
3725   //================================================================================
3726   /*!
3727    * \brief Finds intersection of a curve with a plane
3728    *  \param [in] u1 - parameter of one curve point
3729    *  \param [in] proj1 - projection of the curve point to the plane normal
3730    *  \param [in] u2 - parameter of another curve point
3731    *  \param [in] proj2 - projection of the other curve point to the plane normal
3732    *  \param [in] proj - projection of a point where the curve intersects the plane
3733    *  \param [in] curve - the curve
3734    *  \param [in] axis - the plane normal
3735    *  \param [in] origin - the plane origin
3736    *  \return gp_Pnt - the found intersection point
3737    */
3738   gp_Pnt Hexahedron::findIntPoint( double u1, double proj1,
3739                                    double u2, double proj2,
3740                                    double proj,
3741                                    BRepAdaptor_Curve& curve,
3742                                    const gp_XYZ& axis,
3743                                    const gp_XYZ& origin)
3744   {
3745     double r = (( proj - proj1 ) / ( proj2 - proj1 ));
3746     double u = u1 * ( 1 - r ) + u2 * r;
3747     gp_Pnt p = curve.Value( u );
3748     double newProj =  axis * ( p.XYZ() - origin );
3749     if ( Abs( proj - newProj ) > _grid->_tol / 10. )
3750     {
3751       if ( r > 0.5 )
3752         return findIntPoint( u2, proj2, u, newProj, proj, curve, axis, origin );
3753       else
3754         return findIntPoint( u1, proj2, u, newProj, proj, curve, axis, origin );
3755     }
3756     return p;
3757   }
3758
3759   //================================================================================
3760   /*!
3761    * \brief Returns indices of a hexahedron sub-entities holding a point
3762    *  \param [in] ip - intersection point
3763    *  \param [out] facets - 0-3 facets holding a point
3764    *  \param [out] sub - index of a vertex or an edge holding a point
3765    *  \return int - number of facets holding a point
3766    */
3767   int Hexahedron::getEntity( const E_IntersectPoint* ip, int* facets, int& sub )
3768   {
3769     enum { X = 1, Y = 2, Z = 4 }; // == 001, 010, 100
3770     int nbFacets = 0;
3771     int vertex = 0, edgeMask = 0;
3772
3773     if ( Abs( _grid->_coords[0][ _i   ] - ip->_uvw[0] ) < _grid->_tol ) {
3774       facets[ nbFacets++ ] = SMESH_Block::ID_F0yz;
3775       edgeMask |= X;
3776     }
3777     else if ( Abs( _grid->_coords[0][ _i+1 ] - ip->_uvw[0] ) < _grid->_tol ) {
3778       facets[ nbFacets++ ] = SMESH_Block::ID_F1yz;
3779       vertex   |= X;
3780       edgeMask |= X;
3781     }
3782     if ( Abs( _grid->_coords[1][ _j   ] - ip->_uvw[1] ) < _grid->_tol ) {
3783       facets[ nbFacets++ ] = SMESH_Block::ID_Fx0z;
3784       edgeMask |= Y;
3785     }
3786     else if ( Abs( _grid->_coords[1][ _j+1 ] - ip->_uvw[1] ) < _grid->_tol ) {
3787       facets[ nbFacets++ ] = SMESH_Block::ID_Fx1z;
3788       vertex   |= Y;
3789       edgeMask |= Y;
3790     }
3791     if ( Abs( _grid->_coords[2][ _k   ] - ip->_uvw[2] ) < _grid->_tol ) {
3792       facets[ nbFacets++ ] = SMESH_Block::ID_Fxy0;
3793       edgeMask |= Z;
3794     }
3795     else if ( Abs( _grid->_coords[2][ _k+1 ] - ip->_uvw[2] ) < _grid->_tol ) {
3796       facets[ nbFacets++ ] = SMESH_Block::ID_Fxy1;
3797       vertex   |= Z;
3798       edgeMask |= Z;
3799     }
3800
3801     switch ( nbFacets )
3802     {
3803     case 0: sub = 0;         break;
3804     case 1: sub = facets[0]; break;
3805     case 2: {
3806       const int edge [3][8] = {
3807         { SMESH_Block::ID_E00z, SMESH_Block::ID_E10z,
3808           SMESH_Block::ID_E01z, SMESH_Block::ID_E11z },
3809         { SMESH_Block::ID_E0y0, SMESH_Block::ID_E1y0, 0, 0,
3810           SMESH_Block::ID_E0y1, SMESH_Block::ID_E1y1 },
3811         { SMESH_Block::ID_Ex00, 0, SMESH_Block::ID_Ex10, 0,
3812           SMESH_Block::ID_Ex01, 0, SMESH_Block::ID_Ex11 }
3813       };
3814       switch ( edgeMask ) {
3815       case X | Y: sub = edge[ 0 ][ vertex ]; break;
3816       case X | Z: sub = edge[ 1 ][ vertex ]; break;
3817       default:    sub = edge[ 2 ][ vertex ];
3818       }
3819       break;
3820     }
3821     //case 3:
3822     default:
3823       sub = vertex + SMESH_Block::ID_FirstV;
3824     }
3825
3826     return nbFacets;
3827   }
3828   //================================================================================
3829   /*!
3830    * \brief Adds intersection with an EDGE
3831    */
3832   bool Hexahedron::addIntersection( const E_IntersectPoint* ip,
3833                                     vector< Hexahedron* >&  hexes,
3834                                     int ijk[], int dIJK[] )
3835   {
3836     bool added = false;
3837
3838     size_t hexIndex[4] = {
3839       _grid->CellIndex( ijk[0], ijk[1], ijk[2] ),
3840       dIJK[0] ? _grid->CellIndex( ijk[0]+dIJK[0], ijk[1], ijk[2] ) : -1,
3841       dIJK[1] ? _grid->CellIndex( ijk[0], ijk[1]+dIJK[1], ijk[2] ) : -1,
3842       dIJK[2] ? _grid->CellIndex( ijk[0], ijk[1], ijk[2]+dIJK[2] ) : -1
3843     };
3844     for ( int i = 0; i < 4; ++i )
3845     {
3846       if ( hexIndex[i] < hexes.size() && hexes[ hexIndex[i] ] )
3847       {
3848         Hexahedron* h = hexes[ hexIndex[i] ];
3849         h->_eIntPoints.reserve(2);
3850         h->_eIntPoints.push_back( ip );
3851         added = true;
3852 #ifdef _DEBUG_
3853         // check if ip is really inside the hex
3854         if ( h->isOutParam( ip->_uvw ))
3855           throw SALOME_Exception("ip outside a hex");
3856 #endif
3857       }
3858     }
3859     return added;
3860   }
3861   //================================================================================
3862   /*!
3863    * \brief Finds nodes at a path from one node to another via intersections with EDGEs
3864    */
3865   bool Hexahedron::findChain( _Node*          n1,
3866                               _Node*          n2,
3867                               _Face&          quad,
3868                               vector<_Node*>& chn )
3869   {
3870     chn.clear();
3871     chn.push_back( n1 );
3872     for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
3873       if ( !quad._eIntNodes[ iP ]->IsUsedInFace( &quad ) &&
3874            n1->IsLinked( quad._eIntNodes[ iP ]->_intPoint ) &&
3875            n2->IsLinked( quad._eIntNodes[ iP ]->_intPoint ))
3876       {
3877         chn.push_back( quad._eIntNodes[ iP ]);
3878         chn.push_back( n2 );
3879         quad._eIntNodes[ iP ]->_usedInFace = &quad;
3880         return true;
3881       }
3882     bool found;
3883     do
3884     {
3885       found = false;
3886       for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
3887         if ( !quad._eIntNodes[ iP ]->IsUsedInFace( &quad ) &&
3888              chn.back()->IsLinked( quad._eIntNodes[ iP ]->_intPoint ))
3889         {
3890           chn.push_back( quad._eIntNodes[ iP ]);
3891           found = ( quad._eIntNodes[ iP ]->_usedInFace = &quad );
3892           break;
3893         }
3894     } while ( found && ! chn.back()->IsLinked( n2->_intPoint ) );
3895
3896     if ( chn.back() != n2 && chn.back()->IsLinked( n2->_intPoint ))
3897       chn.push_back( n2 );
3898
3899     return chn.size() > 1;
3900   }
3901   //================================================================================
3902   /*!
3903    * \brief Try to heal a polygon whose ends are not connected
3904    */
3905   bool Hexahedron::closePolygon( _Face* polygon, vector<_Node*>& chainNodes ) const
3906   {
3907     int i = -1, nbLinks = polygon->_links.size();
3908     if ( nbLinks < 3 )
3909       return false;
3910     vector< _OrientedLink > newLinks;
3911     // find a node lying on the same FACE as the last one
3912     _Node*   node = polygon->_links.back().LastNode();
3913     int avoidFace = node->IsLinked( polygon->_links.back().FirstNode()->_intPoint );
3914     for ( i = nbLinks - 2; i >= 0; --i )
3915       if ( node->IsLinked( polygon->_links[i].FirstNode()->_intPoint, avoidFace ))
3916         break;
3917     if ( i >= 0 )
3918     {
3919       for ( ; i < nbLinks; ++i )
3920         newLinks.push_back( polygon->_links[i] );
3921     }
3922     else
3923     {
3924       // find a node lying on the same FACE as the first one
3925       node      = polygon->_links[0].FirstNode();
3926       avoidFace = node->IsLinked( polygon->_links[0].LastNode()->_intPoint );
3927       for ( i = 1; i < nbLinks; ++i )
3928         if ( node->IsLinked( polygon->_links[i].LastNode()->_intPoint, avoidFace ))
3929           break;
3930       if ( i < nbLinks )
3931         for ( nbLinks = i + 1, i = 0; i < nbLinks; ++i )
3932           newLinks.push_back( polygon->_links[i] );
3933     }
3934     if ( newLinks.size() > 1 )
3935     {
3936       polygon->_links.swap( newLinks );
3937       chainNodes.clear();
3938       chainNodes.push_back( polygon->_links.back().LastNode() );
3939       chainNodes.push_back( polygon->_links[0].FirstNode() );
3940       return true;
3941     }
3942     return false;
3943   }
3944   //================================================================================
3945   /*!
3946    * \brief Finds nodes on the same EDGE as the first node of avoidSplit.
3947    *
3948    * This function is for
3949    * 1) a case where an EDGE lies on a quad which lies on a FACE
3950    *    so that a part of quad in ON and another part is IN
3951    * 2) INTERNAL FACE passes through the 1st node of avoidSplit
3952    */
3953   bool Hexahedron::findChainOnEdge( const vector< _OrientedLink >& splits,
3954                                     const _OrientedLink&           prevSplit,
3955                                     const _OrientedLink&           avoidSplit,
3956                                     size_t &                       iS,
3957                                     _Face&                         quad,
3958                                     vector<_Node*>&                chn )
3959   {
3960     _Node* pn1 = prevSplit.FirstNode();
3961     _Node* pn2 = prevSplit.LastNode();
3962     int avoidFace = pn1->IsLinked( pn2->_intPoint ); // FACE under the quad
3963     if ( avoidFace < 1 && pn1->_intPoint )
3964       return false;
3965
3966     _Node* n = 0, *stopNode = avoidSplit.LastNode();
3967
3968     chn.clear();
3969     if ( !quad._eIntNodes.empty() ) // connect pn2 with EDGE intersections
3970     {
3971       chn.push_back( pn2 );
3972       bool found;
3973       do
3974       {
3975         found = false;
3976         for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
3977           if (( !quad._eIntNodes[ iP ]->IsUsedInFace( &quad )) &&
3978               ( chn.back()->IsLinked( quad._eIntNodes[ iP ]->_intPoint, avoidFace )) &&
3979               ( !avoidFace || quad._eIntNodes[ iP ]->IsOnFace( avoidFace )))
3980           {
3981             chn.push_back( quad._eIntNodes[ iP ]);
3982             found = ( quad._eIntNodes[ iP ]->_usedInFace = &quad );
3983             break;
3984           }
3985       } while ( found );
3986       pn2 = chn.back();
3987     }
3988
3989     int i;
3990     for ( i = splits.size()-1; i >= 0; --i ) // connect new pn2 (at _eIntNodes) with a split
3991     {
3992       if ( !splits[i] )
3993         continue;
3994
3995       n = splits[i].LastNode();
3996       if ( n == stopNode )
3997         break;
3998       if (( n != pn1 ) &&
3999           ( n->IsLinked( pn2->_intPoint, avoidFace )) &&
4000           ( !avoidFace || n->IsOnFace( avoidFace )))
4001         break;
4002
4003       n = splits[i].FirstNode();
4004       if ( n == stopNode )
4005         break;
4006       if (( n->IsLinked( pn2->_intPoint, avoidFace )) &&
4007           ( !avoidFace || n->IsOnFace( avoidFace )))
4008         break;
4009       n = 0;
4010     }
4011     if ( n && n != stopNode )
4012     {
4013       if ( chn.empty() )
4014         chn.push_back( pn2 );
4015       chn.push_back( n );
4016       iS = i-1;
4017       return true;
4018     }
4019     else if ( !chn.empty() && chn.back()->_isInternalFlags )
4020     {
4021       // INTERNAL FACE partially cuts the quad
4022       for ( int i = chn.size() - 2; i >= 0; --i )
4023         chn.push_back( chn[ i ]);
4024       return true;
4025     }
4026     return false;
4027   }
4028   //================================================================================
4029   /*!
4030    * \brief Checks transition at the ginen intersection node of a link
4031    */
4032   bool Hexahedron::isOutPoint( _Link& link, int iP,
4033                                SMESH_MesherHelper& helper, const Solid* solid ) const
4034   {
4035     bool isOut = false;
4036
4037     if ( link._fIntNodes[iP]->faces().size() == 1 &&
4038          _grid->IsInternal( link._fIntNodes[iP]->face(0) ))
4039       return false;
4040
4041     const bool moreIntPoints = ( iP+1 < (int) link._fIntNodes.size() );
4042
4043     // get 2 _Node's
4044     _Node* n1 = link._fIntNodes[ iP ];
4045     if ( !n1->Node() )
4046       n1 = link._nodes[0];
4047     _Node* n2 = moreIntPoints ? link._fIntNodes[ iP+1 ] : 0;
4048     if ( !n2 || !n2->Node() )
4049       n2 = link._nodes[1];
4050     if ( !n2->Node() )
4051       return true;
4052
4053     // get all FACEs under n1 and n2
4054     set< TGeomID > faceIDs;
4055     if ( moreIntPoints ) faceIDs.insert( link._fIntNodes[iP+1]->faces().begin(),
4056                                          link._fIntNodes[iP+1]->faces().end() );
4057     if ( n2->_intPoint ) faceIDs.insert( n2->_intPoint->_faceIDs.begin(),
4058                                          n2->_intPoint->_faceIDs.end() );
4059     if ( faceIDs.empty() )
4060       return false; // n2 is inside
4061     if ( n1->_intPoint ) faceIDs.insert( n1->_intPoint->_faceIDs.begin(),
4062                                          n1->_intPoint->_faceIDs.end() );
4063     faceIDs.insert( link._fIntNodes[iP]->faces().begin(),
4064                     link._fIntNodes[iP]->faces().end() );
4065
4066     // get a point between 2 nodes
4067     gp_Pnt p1      = n1->Point();
4068     gp_Pnt p2      = n2->Point();
4069     gp_Pnt pOnLink = 0.8 * p1.XYZ() + 0.2 * p2.XYZ();
4070
4071     TopLoc_Location loc;
4072
4073     set< TGeomID >::iterator faceID = faceIDs.begin();
4074     for ( ; faceID != faceIDs.end(); ++faceID )
4075     {
4076       // project pOnLink on a FACE
4077       if ( *faceID < 1 || !solid->Contains( *faceID )) continue;
4078       const TopoDS_Face& face = TopoDS::Face( _grid->Shape( *faceID ));
4079       GeomAPI_ProjectPointOnSurf& proj = helper.GetProjector( face, loc, 0.1*_grid->_tol );
4080       gp_Pnt testPnt = pOnLink.Transformed( loc.Transformation().Inverted() );
4081       proj.Perform( testPnt );
4082       if ( proj.IsDone() && proj.NbPoints() > 0 )       
4083       {
4084         Standard_Real u,v;
4085         proj.LowerDistanceParameters( u,v );
4086
4087         if ( proj.LowerDistance() <= 0.1 * _grid->_tol )
4088         {
4089           isOut = false;
4090         }
4091         else
4092         {
4093           // find isOut by normals
4094           gp_Dir normal;
4095           if ( GeomLib::NormEstim( BRep_Tool::Surface( face, loc ),
4096                                    gp_Pnt2d( u,v ),
4097                                    0.1*_grid->_tol,
4098                                    normal ) < 3 )
4099           {
4100             if ( solid->Orientation( face ) == TopAbs_REVERSED )
4101               normal.Reverse();
4102             gp_Vec v( proj.NearestPoint(), testPnt );
4103             isOut = ( v * normal > 0 );
4104           }
4105         }
4106         if ( !isOut )
4107         {
4108           // classify a projection
4109           if ( !n1->IsOnFace( *faceID ) || !n2->IsOnFace( *faceID ))
4110           {
4111             BRepTopAdaptor_FClass2d cls( face, Precision::Confusion() );
4112             TopAbs_State state = cls.Perform( gp_Pnt2d( u,v ));
4113             if ( state == TopAbs_OUT )
4114             {
4115               isOut = true;
4116               continue;
4117             }
4118           }
4119           return false;
4120         }
4121       }
4122     }
4123     return isOut;
4124   }
4125   //================================================================================
4126   /*!
4127    * \brief Sort nodes on a FACE
4128    */
4129   void Hexahedron::sortVertexNodes(vector<_Node*>& nodes, _Node* curNode, TGeomID faceID)
4130   {
4131     if ( nodes.size() > 20 ) return;
4132
4133     // get shapes under nodes
4134     TGeomID nShapeIds[20], *nShapeIdsEnd = &nShapeIds[0] + nodes.size();
4135     for ( size_t i = 0; i < nodes.size(); ++i )
4136       if ( !( nShapeIds[i] = nodes[i]->ShapeID() ))
4137         return;
4138
4139     // get shapes of the FACE
4140     const TopoDS_Face&  face = TopoDS::Face( _grid->Shape( faceID ));
4141     list< TopoDS_Edge > edges;
4142     list< int >         nbEdges;
4143     int nbW = SMESH_Block::GetOrderedEdges (face, edges, nbEdges);
4144     if ( nbW > 1 ) {
4145       // select a WIRE - remove EDGEs of irrelevant WIREs from edges
4146       list< TopoDS_Edge >::iterator e = edges.begin(), eEnd = e;
4147       list< int >::iterator nE = nbEdges.begin();
4148       for ( ; nbW > 0; ++nE, --nbW )
4149       {
4150         std::advance( eEnd, *nE );
4151         for ( ; e != eEnd; ++e )
4152           for ( int i = 0; i < 2; ++i )
4153           {
4154             TGeomID id = i==0 ?
4155               _grid->ShapeID( *e ) :
4156               _grid->ShapeID( SMESH_MesherHelper::IthVertex( 0, *e ));
4157             if (( id > 0 ) &&
4158                 ( std::find( &nShapeIds[0], nShapeIdsEnd, id ) != nShapeIdsEnd ))
4159             {
4160               edges.erase( eEnd, edges.end() ); // remove rest wires
4161               e = eEnd = edges.end();
4162               --e;
4163               nbW = 0;
4164               break;
4165             }
4166           }
4167         if ( nbW > 0 )
4168           edges.erase( edges.begin(), eEnd ); // remove a current irrelevant wire
4169       }
4170     }
4171     // rotate edges to have the first one at least partially out of the hexa
4172     list< TopoDS_Edge >::iterator e = edges.begin(), eMidOut = edges.end();
4173     for ( ; e != edges.end(); ++e )
4174     {
4175       if ( !_grid->ShapeID( *e ))
4176         continue;
4177       bool isOut = false;
4178       gp_Pnt p;
4179       double uvw[3], f,l;
4180       for ( int i = 0; i < 2 && !isOut; ++i )
4181       {
4182         if ( i == 0 )
4183         {
4184           TopoDS_Vertex v = SMESH_MesherHelper::IthVertex( 0, *e );
4185           p = BRep_Tool::Pnt( v );
4186         }
4187         else if ( eMidOut == edges.end() )
4188         {
4189           TopLoc_Location loc;
4190           Handle(Geom_Curve) c = BRep_Tool::Curve( *e, loc, f, l);
4191           if ( c.IsNull() ) break;
4192           p = c->Value( 0.5 * ( f + l )).Transformed( loc );
4193         }
4194         else
4195         {
4196           continue;
4197         }
4198
4199         _grid->ComputeUVW( p.XYZ(), uvw );
4200         if ( isOutParam( uvw ))
4201         {
4202           if ( i == 0 )
4203             isOut = true;
4204           else
4205             eMidOut = e;
4206         }
4207       }
4208       if ( isOut )
4209         break;
4210     }
4211     if ( e != edges.end() )
4212       edges.splice( edges.end(), edges, edges.begin(), e );
4213     else if ( eMidOut != edges.end() )
4214       edges.splice( edges.end(), edges, edges.begin(), eMidOut );
4215
4216     // sort nodes according to the order of edges
4217     _Node*  orderNodes   [20];
4218     //TGeomID orderShapeIDs[20];
4219     size_t nbN = 0;
4220     TGeomID id, *pID = 0;
4221     for ( e = edges.begin(); e != edges.end(); ++e )
4222     {
4223       if (( id = _grid->ShapeID( SMESH_MesherHelper::IthVertex( 0, *e ))) &&
4224           (( pID = std::find( &nShapeIds[0], nShapeIdsEnd, id )) != nShapeIdsEnd ))
4225       {
4226         //orderShapeIDs[ nbN ] = id;
4227         orderNodes   [ nbN++ ] = nodes[ pID - &nShapeIds[0] ];
4228         *pID = -1;
4229       }
4230       if (( id = _grid->ShapeID( *e )) &&
4231           (( pID = std::find( &nShapeIds[0], nShapeIdsEnd, id )) != nShapeIdsEnd ))
4232       {
4233         //orderShapeIDs[ nbN ] = id;
4234         orderNodes   [ nbN++ ] = nodes[ pID - &nShapeIds[0] ];
4235         *pID = -1;
4236       }
4237     }
4238     if ( nbN != nodes.size() )
4239       return;
4240
4241     bool reverse = ( orderNodes[0    ]->Point().SquareDistance( curNode->Point() ) >
4242                      orderNodes[nbN-1]->Point().SquareDistance( curNode->Point() ));
4243
4244     for ( size_t i = 0; i < nodes.size(); ++i )
4245       nodes[ i ] = orderNodes[ reverse ? nbN-1-i : i ];
4246   }
4247
4248   //================================================================================
4249   /*!
4250    * \brief Adds computed elements to the mesh
4251    */
4252   int Hexahedron::addVolumes( SMESH_MesherHelper& helper )
4253   {
4254     F_IntersectPoint noIntPnt;
4255     const bool toCheckNodePos = _grid->IsToCheckNodePos();
4256
4257     int nbAdded = 0;
4258     // add elements resulted from hexahedron intersection
4259     for ( _volumeDef* volDef = &_volumeDefs; volDef; volDef = volDef->_next )
4260     {
4261       vector< const SMDS_MeshNode* > nodes( volDef->_nodes.size() );
4262       for ( size_t iN = 0; iN < nodes.size(); ++iN )
4263       {
4264         if ( !( nodes[iN] = volDef->_nodes[iN].Node() ))
4265         {
4266           if ( const E_IntersectPoint* eip = volDef->_nodes[iN].EdgeIntPnt() )
4267           {
4268             nodes[iN] = volDef->_nodes[iN]._intPoint->_node =
4269               helper.AddNode( eip->_point.X(),
4270                               eip->_point.Y(),
4271                               eip->_point.Z() );
4272             if ( _grid->ShapeType( eip->_shapeID ) == TopAbs_VERTEX )
4273               helper.GetMeshDS()->SetNodeOnVertex( nodes[iN], eip->_shapeID );
4274             else
4275               helper.GetMeshDS()->SetNodeOnEdge( nodes[iN], eip->_shapeID );
4276           }
4277           else
4278             throw SALOME_Exception("Bug: no node at intersection point");
4279         }
4280         else if ( volDef->_nodes[iN]._intPoint &&
4281                   volDef->_nodes[iN]._intPoint->_node == volDef->_nodes[iN]._node )
4282         {
4283           // Update position of node at EDGE intersection;
4284           // see comment to _Node::Add( E_IntersectPoint )
4285           SMESHDS_Mesh* mesh = helper.GetMeshDS();
4286           TGeomID    shapeID = volDef->_nodes[iN].EdgeIntPnt()->_shapeID;
4287           mesh->UnSetNodeOnShape( nodes[iN] );
4288           if ( _grid->ShapeType( shapeID ) == TopAbs_VERTEX )
4289             mesh->SetNodeOnVertex( nodes[iN], shapeID );
4290           else
4291             mesh->SetNodeOnEdge( nodes[iN], shapeID );
4292         }
4293         else if ( toCheckNodePos &&
4294                   !nodes[iN]->isMarked() && 
4295                   _grid->ShapeType( nodes[iN]->GetShapeID() ) == TopAbs_FACE )
4296         {
4297           _grid->SetOnShape( nodes[iN], noIntPnt, /*unset=*/true );
4298           nodes[iN]->setIsMarked( true );
4299         }
4300       }
4301
4302       const SMDS_MeshElement* v = 0;
4303       if ( !volDef->_quantities.empty() )
4304       {
4305         v = helper.AddPolyhedralVolume( nodes, volDef->_quantities );
4306       }
4307       else
4308       {
4309         switch ( nodes.size() )
4310         {
4311         case 8: v = helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],
4312                                       nodes[4],nodes[5],nodes[6],nodes[7] );
4313           break;
4314         case 4: v = helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3] );
4315           break;
4316         case 6: v = helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],nodes[4],nodes[5] );
4317           break;
4318         case 5: v = helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],nodes[4] );
4319           break;
4320         }
4321       }
4322       if (( volDef->_volume = v ))
4323       {
4324         helper.GetMeshDS()->SetMeshElementOnShape( v, volDef->_solidID );
4325         ++nbAdded;
4326       }
4327     }
4328
4329     return nbAdded;
4330   }
4331   //================================================================================
4332   /*!
4333    * \brief Return true if the element is in a hole
4334    */
4335   bool Hexahedron::isInHole() const
4336   {
4337     if ( !_vIntNodes.empty() )
4338       return false;
4339
4340     const size_t ijk[3] = { _i, _j, _k };
4341     F_IntersectPoint curIntPnt;
4342
4343     // consider a cell to be in a hole if all links in any direction
4344     // comes OUT of geometry
4345     for ( int iDir = 0; iDir < 3; ++iDir )
4346     {
4347       const vector<double>& coords = _grid->_coords[ iDir ];
4348       LineIndexer               li = _grid->GetLineIndexer( iDir );
4349       li.SetIJK( _i,_j,_k );
4350       size_t lineIndex[4] = { li.LineIndex  (),
4351                               li.LineIndex10(),
4352                               li.LineIndex01(),
4353                               li.LineIndex11() };
4354       bool allLinksOut = true, hasLinks = false;
4355       for ( int iL = 0; iL < 4 && allLinksOut; ++iL ) // loop on 4 links parallel to iDir
4356       {
4357         const _Link& link = _hexLinks[ iL + 4*iDir ];
4358         // check transition of the first node of a link
4359         const F_IntersectPoint* firstIntPnt = 0;
4360         if ( link._nodes[0]->Node() ) // 1st node is a hexa corner
4361         {
4362           curIntPnt._paramOnLine = coords[ ijk[ iDir ]] - coords[0] + _grid->_tol;
4363           const GridLine& line = _grid->_lines[ iDir ][ lineIndex[ iL ]];
4364           multiset< F_IntersectPoint >::const_iterator ip =
4365             line._intPoints.upper_bound( curIntPnt );
4366           --ip;
4367           firstIntPnt = &(*ip);
4368         }
4369         else if ( !link._fIntPoints.empty() )
4370         {
4371           firstIntPnt = link._fIntPoints[0];
4372         }
4373
4374         if ( firstIntPnt )
4375         {
4376           hasLinks = true;
4377           allLinksOut = ( firstIntPnt->_transition == Trans_OUT &&
4378                           !_grid->IsShared( firstIntPnt->_faceIDs[0] ));
4379         }
4380       }
4381       if ( hasLinks && allLinksOut )
4382         return true;
4383     }
4384     return false;
4385   }
4386
4387   //================================================================================
4388   /*!
4389    * \brief Check if a polyherdon has an edge lying on EDGE shared by strange FACE
4390    *        that will be meshed by other algo
4391    */
4392   bool Hexahedron::hasStrangeEdge() const
4393   {
4394     if ( _eIntPoints.size() < 2 )
4395       return false;
4396
4397     TopTools_MapOfShape edges;
4398     for ( size_t i = 0; i < _eIntPoints.size(); ++i )
4399     {
4400       if ( !_grid->IsStrangeEdge( _eIntPoints[i]->_shapeID ))
4401         continue;
4402       const TopoDS_Shape& s = _grid->Shape( _eIntPoints[i]->_shapeID );
4403       if ( s.ShapeType() == TopAbs_EDGE )
4404       {
4405         if ( ! edges.Add( s ))
4406           return true; // an EDGE encounters twice
4407       }
4408       else
4409       {
4410         PShapeIteratorPtr edgeIt = _grid->_helper->GetAncestors( s,
4411                                                                  *_grid->_helper->GetMesh(),
4412                                                                  TopAbs_EDGE );
4413         while ( const TopoDS_Shape* edge = edgeIt->next() )
4414           if ( ! edges.Add( *edge ))
4415             return true; // an EDGE encounters twice
4416       }
4417     }
4418     return false;
4419   }
4420
4421   //================================================================================
4422   /*!
4423    * \brief Return true if a polyhedron passes _sizeThreshold criterion
4424    */
4425   bool Hexahedron::checkPolyhedronSize( bool cutByInternalFace ) const
4426   {
4427     if ( cutByInternalFace && !_grid->_toUseThresholdForInternalFaces )
4428     {
4429       // check if any polygon fully lies on shared/internal FACEs
4430       for ( size_t iP = 0; iP < _polygons.size(); ++iP )
4431       {
4432         const _Face& polygon = _polygons[iP];
4433         if ( polygon._links.empty() )
4434           continue;
4435         bool allNodesInternal = true;
4436         for ( size_t iL = 0; iL < polygon._links.size() &&  allNodesInternal; ++iL )
4437         {
4438           _Node* n = polygon._links[ iL ].FirstNode();
4439           allNodesInternal = (( n->IsCutByInternal() ) ||
4440                               ( n->_intPoint && _grid->IsAnyShared( n->_intPoint->_faceIDs )));
4441         }
4442         if ( allNodesInternal )
4443           return true;
4444       }
4445     }
4446     if ( this->hasStrangeEdge() )
4447       return true;
4448
4449     double volume = 0;
4450     for ( size_t iP = 0; iP < _polygons.size(); ++iP )
4451     {
4452       const _Face& polygon = _polygons[iP];
4453       if ( polygon._links.empty() )
4454         continue;
4455       gp_XYZ area (0,0,0);
4456       gp_XYZ p1 = polygon._links[ 0 ].FirstNode()->Point().XYZ();
4457       for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
4458       {
4459         gp_XYZ p2 = polygon._links[ iL ].LastNode()->Point().XYZ();
4460         area += p1 ^ p2;
4461         p1 = p2;
4462       }
4463       volume += p1 * area;
4464     }
4465     volume /= 6;
4466
4467     double initVolume = _sideLength[0] * _sideLength[1] * _sideLength[2];
4468
4469     return volume > initVolume / _grid->_sizeThreshold;
4470   }
4471   //================================================================================
4472   /*!
4473    * \brief Tries to create a hexahedron
4474    */
4475   bool Hexahedron::addHexa()
4476   {
4477     int nbQuad = 0, iQuad = -1;
4478     for ( size_t i = 0; i < _polygons.size(); ++i )
4479     {
4480       if ( _polygons[i]._links.empty() )
4481         continue;
4482       if ( _polygons[i]._links.size() != 4 )
4483         return false;
4484       ++nbQuad;
4485       if ( iQuad < 0 )
4486         iQuad = i;
4487     }
4488     if ( nbQuad != 6 )
4489       return false;
4490
4491     _Node* nodes[8];
4492     int nbN = 0;
4493     for ( int iL = 0; iL < 4; ++iL )
4494     {
4495       // a base node
4496       nodes[iL] = _polygons[iQuad]._links[iL].FirstNode();
4497       ++nbN;
4498
4499       // find a top node above the base node
4500       _Link* link = _polygons[iQuad]._links[iL]._link;
4501       if ( !link->_faces[0] || !link->_faces[1] )
4502         return debugDumpLink( link );
4503       // a quadrangle sharing <link> with _polygons[iQuad]
4504       _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[iQuad] )];
4505       for ( int i = 0; i < 4; ++i )
4506         if ( quad->_links[i]._link == link )
4507         {
4508           // 1st node of a link opposite to <link> in <quad>
4509           nodes[iL+4] = quad->_links[(i+2)%4].FirstNode();
4510           ++nbN;
4511           break;
4512         }
4513     }
4514     if ( nbN == 8 )
4515       _volumeDefs.Set( &nodes[0], 8 );
4516
4517     return nbN == 8;
4518   }
4519   //================================================================================
4520   /*!
4521    * \brief Tries to create a tetrahedron
4522    */
4523   bool Hexahedron::addTetra()
4524   {
4525     int iTria = -1;
4526     for ( size_t i = 0; i < _polygons.size() && iTria < 0; ++i )
4527       if ( _polygons[i]._links.size() == 3 )
4528         iTria = i;
4529     if ( iTria < 0 )
4530       return false;
4531
4532     _Node* nodes[4];
4533     nodes[0] = _polygons[iTria]._links[0].FirstNode();
4534     nodes[1] = _polygons[iTria]._links[1].FirstNode();
4535     nodes[2] = _polygons[iTria]._links[2].FirstNode();
4536
4537     _Link* link = _polygons[iTria]._links[0]._link;
4538     if ( !link->_faces[0] || !link->_faces[1] )
4539       return debugDumpLink( link );
4540
4541     // a triangle sharing <link> with _polygons[0]
4542     _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[iTria] )];
4543     for ( int i = 0; i < 3; ++i )
4544       if ( tria->_links[i]._link == link )
4545       {
4546         nodes[3] = tria->_links[(i+1)%3].LastNode();
4547         _volumeDefs.Set( &nodes[0], 4 );
4548         return true;
4549       }
4550
4551     return false;
4552   }
4553   //================================================================================
4554   /*!
4555    * \brief Tries to create a pentahedron
4556    */
4557   bool Hexahedron::addPenta()
4558   {
4559     // find a base triangular face
4560     int iTri = -1;
4561     for ( int iF = 0; iF < 5 && iTri < 0; ++iF )
4562       if ( _polygons[ iF ]._links.size() == 3 )
4563         iTri = iF;
4564     if ( iTri < 0 ) return false;
4565
4566     // find nodes
4567     _Node* nodes[6];
4568     int nbN = 0;
4569     for ( int iL = 0; iL < 3; ++iL )
4570     {
4571       // a base node
4572       nodes[iL] = _polygons[ iTri ]._links[iL].FirstNode();
4573       ++nbN;
4574
4575       // find a top node above the base node
4576       _Link* link = _polygons[ iTri ]._links[iL]._link;
4577       if ( !link->_faces[0] || !link->_faces[1] )
4578         return debugDumpLink( link );
4579       // a quadrangle sharing <link> with a base triangle
4580       _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[ iTri ] )];
4581       if ( quad->_links.size() != 4 ) return false;
4582       for ( int i = 0; i < 4; ++i )
4583         if ( quad->_links[i]._link == link )
4584         {
4585           // 1st node of a link opposite to <link> in <quad>
4586           nodes[iL+3] = quad->_links[(i+2)%4].FirstNode();
4587           ++nbN;
4588           break;
4589         }
4590     }
4591     if ( nbN == 6 )
4592       _volumeDefs.Set( &nodes[0], 6 );
4593
4594     return ( nbN == 6 );
4595   }
4596   //================================================================================
4597   /*!
4598    * \brief Tries to create a pyramid
4599    */
4600   bool Hexahedron::addPyra()
4601   {
4602     // find a base quadrangle
4603     int iQuad = -1;
4604     for ( int iF = 0; iF < 5 && iQuad < 0; ++iF )
4605       if ( _polygons[ iF ]._links.size() == 4 )
4606         iQuad = iF;
4607     if ( iQuad < 0 ) return false;
4608
4609     // find nodes
4610     _Node* nodes[5];
4611     nodes[0] = _polygons[iQuad]._links[0].FirstNode();
4612     nodes[1] = _polygons[iQuad]._links[1].FirstNode();
4613     nodes[2] = _polygons[iQuad]._links[2].FirstNode();
4614     nodes[3] = _polygons[iQuad]._links[3].FirstNode();
4615
4616     _Link* link = _polygons[iQuad]._links[0]._link;
4617     if ( !link->_faces[0] || !link->_faces[1] )
4618       return debugDumpLink( link );
4619
4620     // a triangle sharing <link> with a base quadrangle
4621     _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[ iQuad ] )];
4622     if ( tria->_links.size() != 3 ) return false;
4623     for ( int i = 0; i < 3; ++i )
4624       if ( tria->_links[i]._link == link )
4625       {
4626         nodes[4] = tria->_links[(i+1)%3].LastNode();
4627         _volumeDefs.Set( &nodes[0], 5 );
4628         return true;
4629       }
4630
4631     return false;
4632   }
4633   //================================================================================
4634   /*!
4635    * \brief Dump a link and return \c false
4636    */
4637   bool Hexahedron::debugDumpLink( Hexahedron::_Link* link )
4638   {
4639 #ifdef _DEBUG_
4640     gp_Pnt p1 = link->_nodes[0]->Point(), p2 = link->_nodes[1]->Point();
4641     cout << "BUG: not shared link. IKJ = ( "<< _i << " " << _j << " " << _k << " )" << endl
4642          << "n1 (" << p1.X() << ", "<< p1.Y() << ", "<< p1.Z() << " )" << endl
4643          << "n2 (" << p2.X() << ", "<< p2.Y() << ", "<< p2.Z() << " )" << endl;
4644 #endif
4645     return false;
4646   }
4647   //================================================================================
4648   /*!
4649    * \brief Classify a point by grid parameters
4650    */
4651   bool Hexahedron::isOutParam(const double uvw[3]) const
4652   {
4653     return (( _grid->_coords[0][ _i   ] - _grid->_tol > uvw[0] ) ||
4654             ( _grid->_coords[0][ _i+1 ] + _grid->_tol < uvw[0] ) ||
4655             ( _grid->_coords[1][ _j   ] - _grid->_tol > uvw[1] ) ||
4656             ( _grid->_coords[1][ _j+1 ] + _grid->_tol < uvw[1] ) ||
4657             ( _grid->_coords[2][ _k   ] - _grid->_tol > uvw[2] ) ||
4658             ( _grid->_coords[2][ _k+1 ] + _grid->_tol < uvw[2] ));
4659   }
4660   //================================================================================
4661   /*!
4662    * \brief Divide a polygon into triangles and modify accordingly an adjacent polyhedron
4663    */
4664   void splitPolygon( const SMDS_MeshElement*         polygon,
4665                      SMDS_VolumeTool &               volume,
4666                      const int                       facetIndex,
4667                      const TGeomID                   faceID,
4668                      const TGeomID                   solidID,
4669                      SMESH_MeshEditor::ElemFeatures& face,
4670                      SMESH_MeshEditor&               editor,
4671                      const bool                      reinitVolume)
4672   {
4673     SMESH_MeshAlgos::Triangulate divider(/*optimize=*/false);
4674     int nbTrias = divider.GetTriangles( polygon, face.myNodes );
4675     face.myNodes.resize( nbTrias * 3 );
4676
4677     SMESH_MeshEditor::ElemFeatures newVolumeDef;
4678     newVolumeDef.Init( volume.Element() );
4679     newVolumeDef.SetID( volume.Element()->GetID() );
4680
4681     newVolumeDef.myPolyhedQuantities.reserve( volume.NbFaces() + nbTrias );
4682     newVolumeDef.myNodes.reserve( volume.NbNodes() + nbTrias * 3 );
4683
4684     SMESHDS_Mesh* meshDS = editor.GetMeshDS();
4685     SMDS_MeshElement* newTriangle;
4686     for ( int iF = 0, nF = volume.NbFaces(); iF < nF; iF++ )
4687     {
4688       if ( iF == facetIndex )
4689       {
4690         newVolumeDef.myPolyhedQuantities.push_back( 3 );
4691         newVolumeDef.myNodes.insert( newVolumeDef.myNodes.end(),
4692                                      face.myNodes.begin(),
4693                                      face.myNodes.begin() + 3 );
4694         meshDS->RemoveFreeElement( polygon, 0, false );
4695         newTriangle = meshDS->AddFace( face.myNodes[0], face.myNodes[1], face.myNodes[2] );
4696         meshDS->SetMeshElementOnShape( newTriangle, faceID );
4697       }
4698       else
4699       {
4700         const SMDS_MeshNode** nn = volume.GetFaceNodes( iF );
4701         const size_t nbFaceNodes = volume.NbFaceNodes ( iF );
4702         newVolumeDef.myPolyhedQuantities.push_back( nbFaceNodes );
4703         newVolumeDef.myNodes.insert( newVolumeDef.myNodes.end(), nn, nn + nbFaceNodes );
4704       }
4705     }
4706
4707     for ( size_t iN = 3; iN < face.myNodes.size(); iN += 3 )
4708     {
4709       newVolumeDef.myPolyhedQuantities.push_back( 3 );
4710       newVolumeDef.myNodes.insert( newVolumeDef.myNodes.end(),
4711                                    face.myNodes.begin() + iN,
4712                                    face.myNodes.begin() + iN + 3 );
4713       newTriangle = meshDS->AddFace( face.myNodes[iN], face.myNodes[iN+1], face.myNodes[iN+2] );
4714       meshDS->SetMeshElementOnShape( newTriangle, faceID );
4715     }
4716
4717     meshDS->RemoveFreeElement( volume.Element(), 0, false );
4718     SMDS_MeshElement* newVolume = editor.AddElement( newVolumeDef.myNodes, newVolumeDef );
4719     meshDS->SetMeshElementOnShape( newVolume, solidID );
4720
4721     if ( reinitVolume )
4722     {
4723       volume.Set( 0 );
4724       volume.Set( newVolume );
4725     }
4726     return;
4727   }
4728   //================================================================================
4729   /*!
4730    * \brief Create mesh faces at free facets
4731    */
4732   void Hexahedron::addFaces( SMESH_MesherHelper&                       helper,
4733                              const vector< const SMDS_MeshElement* > & boundaryVolumes )
4734   {
4735     if ( !_grid->_toCreateFaces )
4736       return;
4737
4738     SMDS_VolumeTool vTool;
4739     vector<int> bndFacets;
4740     SMESH_MeshEditor editor( helper.GetMesh() );
4741     SMESH_MeshEditor::ElemFeatures face( SMDSAbs_Face );
4742     SMESHDS_Mesh* meshDS = helper.GetMeshDS();
4743
4744     // check if there are internal or shared FACEs
4745     bool hasInternal = ( !_grid->_geometry.IsOneSolid() ||
4746                          _grid->_geometry._soleSolid.HasInternalFaces() );
4747
4748     for ( size_t iV = 0; iV < boundaryVolumes.size(); ++iV )
4749     {
4750       if ( !vTool.Set( boundaryVolumes[ iV ]))
4751         continue;
4752
4753       TGeomID solidID = vTool.Element()->GetShapeID();
4754       Solid *   solid = _grid->GetOneOfSolids( solidID );
4755
4756       // find boundary facets
4757
4758       bndFacets.clear();
4759       for ( int iF = 0, n = vTool.NbFaces(); iF < n; iF++ )
4760       {
4761         bool isBoundary = vTool.IsFreeFace( iF );
4762         if ( isBoundary )
4763         {
4764           bndFacets.push_back( iF );
4765         }
4766         else if ( hasInternal )
4767         {
4768           // check if all nodes are on internal/shared FACEs
4769           isBoundary = true;
4770           const SMDS_MeshNode** nn = vTool.GetFaceNodes( iF );
4771           const size_t nbFaceNodes = vTool.NbFaceNodes ( iF );
4772           for ( size_t iN = 0; iN < nbFaceNodes &&  isBoundary; ++iN )
4773             isBoundary = ( nn[ iN ]->GetShapeID() != solidID );
4774           if ( isBoundary )
4775             bndFacets.push_back( -( iF+1 )); // !!! minus ==> to check the FACE
4776         }
4777       }
4778       if ( bndFacets.empty() )
4779         continue;
4780
4781       // create faces
4782
4783       if ( !vTool.IsPoly() )
4784         vTool.SetExternalNormal();
4785       for ( size_t i = 0; i < bndFacets.size(); ++i ) // loop on boundary facets
4786       {
4787         const bool    isBoundary = ( bndFacets[i] >= 0 );
4788         const int         iFacet = isBoundary ? bndFacets[i] : -bndFacets[i]-1;
4789         const SMDS_MeshNode** nn = vTool.GetFaceNodes( iFacet );
4790         const size_t nbFaceNodes = vTool.NbFaceNodes ( iFacet );
4791         face.myNodes.assign( nn, nn + nbFaceNodes );
4792
4793         TGeomID faceID = 0;
4794         const SMDS_MeshElement* existFace = 0, *newFace = 0;
4795
4796         if (( existFace = meshDS->FindElement( face.myNodes, SMDSAbs_Face )))
4797         {
4798           if ( existFace->isMarked() )
4799             continue; // created by this method
4800           faceID = existFace->GetShapeID();
4801         }
4802         else
4803         {
4804           // look for a supporting FACE
4805           for ( size_t iN = 0; iN < nbFaceNodes &&  !faceID; ++iN ) // look for a node on FACE
4806           {
4807             if ( nn[ iN ]->GetPosition()->GetDim() == 2 )
4808               faceID = nn[ iN ]->GetShapeID();
4809           }
4810           for ( size_t iN = 0; iN < nbFaceNodes &&  !faceID; ++iN )
4811           {
4812             // look for a father FACE of EDGEs and VERTEXes
4813             const TopoDS_Shape& s1 = _grid->Shape( nn[ iN   ]->GetShapeID() );
4814             const TopoDS_Shape& s2 = _grid->Shape( nn[ iN+1 ]->GetShapeID() );
4815             if ( s1 != s2 && s1.ShapeType() == TopAbs_EDGE && s2.ShapeType() == TopAbs_EDGE )
4816             {
4817               TopoDS_Shape f = helper.GetCommonAncestor( s1, s2, *helper.GetMesh(), TopAbs_FACE );
4818               if ( !f.IsNull() )
4819                 faceID = _grid->ShapeID( f );
4820             }
4821           }
4822
4823           bool toCheckFace = faceID && (( !isBoundary ) ||
4824                                         ( hasInternal && _grid->_toUseThresholdForInternalFaces ));
4825           if ( toCheckFace ) // check if all nodes are on the found FACE
4826           {
4827             SMESH_subMesh* faceSM = helper.GetMesh()->GetSubMeshContaining( faceID );
4828             for ( size_t iN = 0; iN < nbFaceNodes &&  faceID; ++iN )
4829             {
4830               TGeomID subID = nn[ iN ]->GetShapeID();
4831               if ( subID != faceID && !faceSM->DependsOn( subID ))
4832                 faceID = 0;
4833             }
4834             if ( !faceID && !isBoundary )
4835               continue;
4836           }
4837         }
4838         // orient a new face according to supporting FACE orientation in shape_to_mesh
4839         if ( !solid->IsOutsideOriented( faceID ))
4840         {
4841           if ( existFace )
4842             editor.Reorient( existFace );
4843           else
4844             std::reverse( face.myNodes.begin(), face.myNodes.end() );
4845         }
4846
4847         if ( ! ( newFace = existFace ))
4848         {
4849           face.SetPoly( nbFaceNodes > 4 );
4850           newFace = editor.AddElement( face.myNodes, face );
4851           if ( !newFace )
4852             continue;
4853           newFace->setIsMarked( true ); // to distinguish from face created in getBoundaryElems()
4854         }
4855
4856         if ( faceID && _grid->IsBoundaryFace( faceID )) // face is not shared
4857         {
4858           // set newFace to the found FACE provided that it fully lies on the FACE
4859           for ( size_t iN = 0; iN < nbFaceNodes &&  faceID; ++iN )
4860             if ( nn[iN]->GetShapeID() == solidID )
4861             {
4862               if ( existFace )
4863                 meshDS->UnSetMeshElementOnShape( existFace, _grid->Shape( faceID ));
4864               faceID = 0;
4865             }
4866         }
4867
4868         // split a polygon that will be used by other 3D algorithm
4869         if ( faceID && nbFaceNodes > 4 &&
4870              !_grid->IsInternal( faceID ) &&
4871              !_grid->IsShared( faceID ) &&
4872              !_grid->IsBoundaryFace( faceID ))
4873         {
4874           splitPolygon( newFace, vTool, iFacet, faceID, solidID,
4875                         face, editor, i+1 < bndFacets.size() );
4876         }
4877         else
4878         {
4879           if ( faceID )
4880             meshDS->SetMeshElementOnShape( newFace, faceID );
4881           else
4882             meshDS->SetMeshElementOnShape( newFace, solidID );
4883         }
4884       } // loop on bndFacets
4885     } // loop on boundaryVolumes
4886
4887
4888     // Orient coherently mesh faces on INTERNAL FACEs
4889
4890     if ( hasInternal )
4891     {
4892       TopExp_Explorer faceExp( _grid->_geometry._mainShape, TopAbs_FACE );
4893       for ( ; faceExp.More(); faceExp.Next() )
4894       {
4895         if ( faceExp.Current().Orientation() != TopAbs_INTERNAL )
4896           continue;
4897
4898         SMESHDS_SubMesh* sm = meshDS->MeshElements( faceExp.Current() );
4899         if ( !sm ) continue;
4900
4901         TIDSortedElemSet facesToOrient;
4902         for ( SMDS_ElemIteratorPtr fIt = sm->GetElements(); fIt->more(); )
4903           facesToOrient.insert( facesToOrient.end(), fIt->next() );
4904         if ( facesToOrient.size() < 2 )
4905           continue;
4906
4907         gp_Dir direction(1,0,0);
4908         const SMDS_MeshElement* anyFace = *facesToOrient.begin();
4909         editor.Reorient2D( facesToOrient, direction, anyFace );
4910       }
4911     }
4912     return;
4913   }
4914
4915   //================================================================================
4916   /*!
4917    * \brief Create mesh segments.
4918    */
4919   void Hexahedron::addSegments( SMESH_MesherHelper&                      helper,
4920                                 const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap )
4921   {
4922     SMESHDS_Mesh* mesh = helper.GetMeshDS();
4923
4924     std::vector<const SMDS_MeshNode*> nodes;
4925     std::vector<const SMDS_MeshElement *> elems;
4926     map< TGeomID, vector< TGeomID > >::const_iterator e2ff = edge2faceIDsMap.begin();
4927     for ( ; e2ff != edge2faceIDsMap.end(); ++e2ff )
4928     {
4929       const TopoDS_Edge& edge = TopoDS::Edge( _grid->Shape( e2ff->first ));
4930       const TopoDS_Face& face = TopoDS::Face( _grid->Shape( e2ff->second[0] ));
4931       StdMeshers_FaceSide side( face, edge, helper.GetMesh(), /*isFwd=*/true, /*skipMed=*/true );
4932       nodes = side.GetOrderedNodes();
4933
4934       elems.clear();
4935       if ( nodes.size() == 2 )
4936         // check that there is an element connecting two nodes
4937         if ( !mesh->GetElementsByNodes( nodes, elems ))
4938           continue;
4939
4940       for ( size_t i = 1; i < nodes.size(); i++ )
4941       {
4942         SMDS_MeshElement* segment = mesh->AddEdge( nodes[i-1], nodes[i] );
4943         mesh->SetMeshElementOnShape( segment, e2ff->first );
4944       }
4945     }
4946     return;
4947   }
4948
4949   //================================================================================
4950   /*!
4951    * \brief Return created volumes and volumes that can have free facet because of
4952    *        skipped small volume. Also create mesh faces on free facets
4953    *        of adjacent not-cut volumes id the result volume is too small.
4954    */
4955   void Hexahedron::getBoundaryElems( vector< const SMDS_MeshElement* > & boundaryElems )
4956   {
4957     if ( _hasTooSmall /*|| _volumeDefs.IsEmpty()*/ )
4958     {
4959       // create faces around a missing small volume
4960       TGeomID faceID = 0;
4961       SMESH_MeshEditor editor( _grid->_helper->GetMesh() );
4962       SMESH_MeshEditor::ElemFeatures polygon( SMDSAbs_Face );
4963       SMESHDS_Mesh* meshDS = _grid->_helper->GetMeshDS();
4964       std::vector<const SMDS_MeshElement *> adjVolumes(2);
4965       for ( size_t iF = 0; iF < _polygons.size(); ++iF )
4966       {
4967         const size_t nbLinks = _polygons[ iF ]._links.size();
4968         if ( nbLinks != 4 ) continue;
4969         polygon.myNodes.resize( nbLinks );
4970         polygon.myNodes.back() = 0;
4971         for ( size_t iL = 0, iN = nbLinks - 1; iL < nbLinks; ++iL, --iN )
4972           if ( ! ( polygon.myNodes[iN] = _polygons[ iF ]._links[ iL ].FirstNode()->Node() ))
4973             break;
4974         if ( !polygon.myNodes.back() )
4975           continue;
4976
4977         meshDS->GetElementsByNodes( polygon.myNodes, adjVolumes, SMDSAbs_Volume );
4978         if ( adjVolumes.size() != 1 )
4979           continue;
4980         if ( !adjVolumes[0]->isMarked() )
4981         {
4982           boundaryElems.push_back( adjVolumes[0] );
4983           adjVolumes[0]->setIsMarked( true );
4984         }
4985
4986         bool sameShape = true;
4987         TGeomID shapeID = polygon.myNodes[0]->GetShapeID();
4988         for ( size_t i = 1; i < polygon.myNodes.size() && sameShape; ++i )
4989           sameShape = ( shapeID == polygon.myNodes[i]->GetShapeID() );
4990
4991         if ( !sameShape || !_grid->IsSolid( shapeID ))
4992           continue; // some of shapes must be FACE
4993
4994         if ( !faceID )
4995         {
4996           faceID = getAnyFace();
4997           if ( !faceID )
4998             break;
4999           if ( _grid->IsInternal( faceID ) ||
5000                _grid->IsShared( faceID ) ||
5001                _grid->IsBoundaryFace( faceID ))
5002             break; // create only if a new face will be used by other 3D algo
5003         }
5004
5005         Solid * solid = _grid->GetOneOfSolids( adjVolumes[0]->GetShapeID() );
5006         if ( !solid->IsOutsideOriented( faceID ))
5007           std::reverse( polygon.myNodes.begin(), polygon.myNodes.end() );
5008
5009         //polygon.SetPoly( polygon.myNodes.size() > 4 );
5010         const SMDS_MeshElement* newFace = editor.AddElement( polygon.myNodes, polygon );
5011         meshDS->SetMeshElementOnShape( newFace, faceID );
5012       }
5013     }
5014
5015     // return created volumes
5016     for ( _volumeDef* volDef = &_volumeDefs; volDef; volDef = volDef->_next )
5017     {
5018       if ( volDef->_volume && !volDef->_volume->isMarked() )
5019       {
5020         volDef->_volume->setIsMarked( true );
5021         boundaryElems.push_back( volDef->_volume );
5022
5023         if ( _grid->IsToCheckNodePos() ) // un-mark nodes marked in addVolumes()
5024           for ( size_t iN = 0; iN < volDef->_nodes.size(); ++iN )
5025             volDef->_nodes[iN].Node()->setIsMarked( false );
5026       }
5027     }
5028   }
5029
5030   //================================================================================
5031   /*!
5032    * \brief Set to _hexLinks a next portion of splits located on one side of INTERNAL FACEs
5033    */
5034   bool Hexahedron::_SplitIterator::Next()
5035   {
5036     if ( _iterationNb > 0 )
5037       // count used splits
5038       for ( size_t i = 0; i < _splits.size(); ++i )
5039       {
5040         if ( _splits[i]._iCheckIteration == _iterationNb )
5041         {
5042           _splits[i]._isUsed = _splits[i]._checkedSplit->_faces[1];
5043           _nbUsed += _splits[i]._isUsed;
5044         }
5045         if ( !More() )
5046           return false;
5047       }
5048
5049     ++_iterationNb;
5050
5051     bool toTestUsed = ( _nbChecked >= _splits.size() );
5052     if ( toTestUsed )
5053     {
5054       // all splits are checked; find all not used splits
5055       for ( size_t i = 0; i < _splits.size(); ++i )
5056         if ( !_splits[i].IsCheckedOrUsed( toTestUsed ))
5057           _splits[i]._iCheckIteration = _iterationNb;
5058
5059       _nbUsed = _splits.size(); // to stop iteration
5060     }
5061     else
5062     {
5063       // get any not used/checked split to start from
5064       _freeNodes.clear();
5065       for ( size_t i = 0; i < _splits.size(); ++i )
5066       {
5067         if ( !_splits[i].IsCheckedOrUsed( toTestUsed ))
5068         {
5069           _freeNodes.push_back( _splits[i]._nodes[0] );
5070           _freeNodes.push_back( _splits[i]._nodes[1] );
5071           _splits[i]._iCheckIteration = _iterationNb;
5072           break;
5073         }
5074       }
5075       // find splits connected to the start one via _freeNodes
5076       for ( size_t iN = 0; iN < _freeNodes.size(); ++iN )
5077       {
5078         for ( size_t iS = 0; iS < _splits.size(); ++iS )
5079         {
5080           if ( _splits[iS].IsCheckedOrUsed( toTestUsed ))
5081             continue;
5082           int iN2 = -1;
5083           if (      _freeNodes[iN] == _splits[iS]._nodes[0] )
5084             iN2 = 1;
5085           else if ( _freeNodes[iN] == _splits[iS]._nodes[1] )
5086             iN2 = 0;
5087           else
5088             continue;
5089           if ( _freeNodes[iN]->_isInternalFlags > 0 )
5090           {
5091             if ( _splits[iS]._nodes[ iN2 ]->_isInternalFlags == 0 )
5092               continue;
5093             if ( !_splits[iS]._nodes[ iN2 ]->IsLinked( _freeNodes[iN]->_intPoint ))
5094               continue;
5095           }
5096           _splits[iS]._iCheckIteration = _iterationNb;
5097           _freeNodes.push_back( _splits[iS]._nodes[ iN2 ]);
5098         }
5099       }
5100     }
5101     // set splits to hex links
5102
5103     for ( int iL = 0; iL < 12; ++iL )
5104       _hexLinks[ iL ]._splits.clear();
5105
5106     _Link split;
5107     for ( size_t i = 0; i < _splits.size(); ++i )
5108     {
5109       if ( _splits[i]._iCheckIteration == _iterationNb )
5110       {
5111         split._nodes[0] = _splits[i]._nodes[0];
5112         split._nodes[1] = _splits[i]._nodes[1];
5113         _Link & hexLink = _hexLinks[ _splits[i]._linkID ];
5114         hexLink._splits.push_back( split );
5115         _splits[i]._checkedSplit = & hexLink._splits.back();
5116         ++_nbChecked;
5117       }
5118     }
5119     return More();
5120   }
5121
5122   //================================================================================
5123   /*!
5124    * \brief computes exact bounding box with axes parallel to given ones
5125    */
5126   //================================================================================
5127
5128   void getExactBndBox( const vector< TopoDS_Shape >& faceVec,
5129                        const double*                 axesDirs,
5130                        Bnd_Box&                      shapeBox )
5131   {
5132     BRep_Builder b;
5133     TopoDS_Compound allFacesComp;
5134     b.MakeCompound( allFacesComp );
5135     for ( size_t iF = 0; iF < faceVec.size(); ++iF )
5136       b.Add( allFacesComp, faceVec[ iF ] );
5137
5138     double sP[6]; // aXmin, aYmin, aZmin, aXmax, aYmax, aZmax
5139     shapeBox.Get(sP[0],sP[1],sP[2],sP[3],sP[4],sP[5]);
5140     double farDist = 0;
5141     for ( int i = 0; i < 6; ++i )
5142       farDist = Max( farDist, 10 * sP[i] );
5143
5144     gp_XYZ axis[3] = { gp_XYZ( axesDirs[0], axesDirs[1], axesDirs[2] ),
5145                        gp_XYZ( axesDirs[3], axesDirs[4], axesDirs[5] ),
5146                        gp_XYZ( axesDirs[6], axesDirs[7], axesDirs[8] ) };
5147     axis[0].Normalize();
5148     axis[1].Normalize();
5149     axis[2].Normalize();
5150
5151     gp_Mat basis( axis[0], axis[1], axis[2] );
5152     gp_Mat bi = basis.Inverted();
5153
5154     gp_Pnt pMin, pMax;
5155     for ( int iDir = 0; iDir < 3; ++iDir )
5156     {
5157       gp_XYZ axis0 = axis[ iDir ];
5158       gp_XYZ axis1 = axis[ ( iDir + 1 ) % 3 ];
5159       gp_XYZ axis2 = axis[ ( iDir + 2 ) % 3 ];
5160       for ( int isMax = 0; isMax < 2; ++isMax )
5161       {
5162         double shift = isMax ? farDist : -farDist;
5163         gp_XYZ orig = shift * axis0;
5164         gp_XYZ norm = axis1 ^ axis2;
5165         gp_Pln pln( orig, norm );
5166         norm = pln.Axis().Direction().XYZ();
5167         BRepBuilderAPI_MakeFace plane( pln, -farDist, farDist, -farDist, farDist );
5168
5169         gp_Pnt& pAxis = isMax ? pMax : pMin;
5170         gp_Pnt pPlane, pFaces;
5171         double dist = GEOMUtils::GetMinDistance( plane, allFacesComp, pPlane, pFaces );
5172         if ( dist < 0 )
5173         {
5174           Bnd_B3d bb;
5175           gp_XYZ corner;
5176           for ( int i = 0; i < 2; ++i ) {
5177             corner.SetCoord( 1, sP[ i*3 ]);
5178             for ( int j = 0; j < 2; ++j ) {
5179               corner.SetCoord( 2, sP[ i*3 + 1 ]);
5180               for ( int k = 0; k < 2; ++k )
5181               {
5182                 corner.SetCoord( 3, sP[ i*3 + 2 ]);
5183                 corner *= bi;
5184                 bb.Add( corner );
5185               }
5186             }
5187           }
5188           corner = isMax ? bb.CornerMax() : bb.CornerMin();
5189           pAxis.SetCoord( iDir+1, corner.Coord( iDir+1 ));
5190         }
5191         else
5192         {
5193           gp_XYZ pf = pFaces.XYZ() * bi;
5194           pAxis.SetCoord( iDir+1, pf.Coord( iDir+1 ) );
5195         }
5196       }
5197     } // loop on 3 axes
5198
5199     shapeBox.SetVoid();
5200     shapeBox.Add( pMin );
5201     shapeBox.Add( pMax );
5202
5203     return;
5204   }
5205
5206 } // namespace
5207
5208 //=============================================================================
5209 /*!
5210  * \brief Generates 3D structured Cartesian mesh in the internal part of
5211  * solid shapes and polyhedral volumes near the shape boundary.
5212  *  \param theMesh - mesh to fill in
5213  *  \param theShape - a compound of all SOLIDs to mesh
5214  *  \retval bool - true in case of success
5215  */
5216 //=============================================================================
5217
5218 bool StdMeshers_Cartesian_3D::Compute(SMESH_Mesh &         theMesh,
5219                                       const TopoDS_Shape & theShape)
5220 {
5221   // The algorithm generates the mesh in following steps:
5222
5223   // 1) Intersection of grid lines with the geometry boundary.
5224   // This step allows to find out if a given node of the initial grid is
5225   // inside or outside the geometry.
5226
5227   // 2) For each cell of the grid, check how many of it's nodes are outside
5228   // of the geometry boundary. Depending on a result of this check
5229   // - skip a cell, if all it's nodes are outside
5230   // - skip a cell, if it is too small according to the size threshold
5231   // - add a hexahedron in the mesh, if all nodes are inside
5232   // - add a polyhedron in the mesh, if some nodes are inside and some outside
5233
5234   _computeCanceled = false;
5235
5236   SMESH_MesherHelper helper( theMesh );
5237   SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
5238
5239   try
5240   {
5241     Grid grid;
5242     grid._helper                         = &helper;
5243     grid._toAddEdges                     = _hyp->GetToAddEdges();
5244     grid._toCreateFaces                  = _hyp->GetToCreateFaces();
5245     grid._toConsiderInternalFaces        = _hyp->GetToConsiderInternalFaces();
5246     grid._toUseThresholdForInternalFaces = _hyp->GetToUseThresholdForInternalFaces();
5247     grid._sizeThreshold                  = _hyp->GetSizeThreshold();
5248     grid.InitGeometry( theShape );
5249
5250     vector< TopoDS_Shape > faceVec;
5251     {
5252       TopTools_MapOfShape faceMap;
5253       TopExp_Explorer fExp;
5254       for ( fExp.Init( theShape, TopAbs_FACE ); fExp.More(); fExp.Next() )
5255       {
5256         bool isNewFace = faceMap.Add( fExp.Current() );
5257         if ( !grid._toConsiderInternalFaces )
5258           if ( !isNewFace || fExp.Current().Orientation() == TopAbs_INTERNAL )
5259             // remove an internal face
5260             faceMap.Remove( fExp.Current() );
5261       }
5262       faceVec.reserve( faceMap.Extent() );
5263       faceVec.assign( faceMap.cbegin(), faceMap.cend() );
5264     }
5265     vector<FaceGridIntersector> facesItersectors( faceVec.size() );
5266     Bnd_Box shapeBox;
5267     for ( size_t i = 0; i < faceVec.size(); ++i )
5268     {
5269       facesItersectors[i]._face   = TopoDS::Face( faceVec[i] );
5270       facesItersectors[i]._faceID = grid.ShapeID( faceVec[i] );
5271       facesItersectors[i]._grid   = &grid;
5272       shapeBox.Add( facesItersectors[i].GetFaceBndBox() );
5273     }
5274     getExactBndBox( faceVec, _hyp->GetAxisDirs(), shapeBox );
5275
5276
5277     vector<double> xCoords, yCoords, zCoords;
5278     _hyp->GetCoordinates( xCoords, yCoords, zCoords, shapeBox );
5279
5280     grid.SetCoordinates( xCoords, yCoords, zCoords, _hyp->GetAxisDirs(), shapeBox );
5281
5282     if ( _computeCanceled ) return false;
5283
5284 #ifdef WITH_TBB
5285     { // copy partner faces and curves of not thread-safe types
5286       set< const Standard_Transient* > tshapes;
5287       BRepBuilderAPI_Copy copier;
5288       for ( size_t i = 0; i < facesItersectors.size(); ++i )
5289       {
5290         if ( !facesItersectors[i].IsThreadSafe( tshapes ))
5291         {
5292           copier.Perform( facesItersectors[i]._face );
5293           facesItersectors[i]._face = TopoDS::Face( copier );
5294         }
5295       }
5296     }
5297     // Intersection of grid lines with the geometry boundary.
5298     tbb::parallel_for ( tbb::blocked_range<size_t>( 0, facesItersectors.size() ),
5299                         ParallelIntersector( facesItersectors ),
5300                         tbb::simple_partitioner());
5301 #else
5302     for ( size_t i = 0; i < facesItersectors.size(); ++i )
5303       facesItersectors[i].Intersect();
5304 #endif
5305
5306     // put intersection points onto the GridLine's; this is done after intersection
5307     // to avoid contention of facesItersectors for writing into the same GridLine
5308     // in case of parallel work of facesItersectors
5309     for ( size_t i = 0; i < facesItersectors.size(); ++i )
5310       facesItersectors[i].StoreIntersections();
5311
5312     if ( _computeCanceled ) return false;
5313
5314     // create nodes on the geometry
5315     grid.ComputeNodes( helper );
5316
5317     if ( _computeCanceled ) return false;
5318
5319     // get EDGEs to take into account
5320     map< TGeomID, vector< TGeomID > > edge2faceIDsMap;
5321     grid.GetEdgesToImplement( edge2faceIDsMap, theShape, faceVec );
5322
5323     // create volume elements
5324     Hexahedron hex( &grid );
5325     int nbAdded = hex.MakeElements( helper, edge2faceIDsMap );
5326
5327     if ( nbAdded > 0 )
5328     {
5329       if ( !grid._toConsiderInternalFaces )
5330       {
5331         // make all SOLIDs computed
5332         TopExp_Explorer solidExp( theShape, TopAbs_SOLID );
5333         if ( SMESHDS_SubMesh* sm1 = meshDS->MeshElements( solidExp.Current()) )
5334         {
5335           SMDS_ElemIteratorPtr volIt = sm1->GetElements();
5336           for ( ; solidExp.More() && volIt->more(); solidExp.Next() )
5337           {
5338             const SMDS_MeshElement* vol = volIt->next();
5339             sm1->RemoveElement( vol );
5340             meshDS->SetMeshElementOnShape( vol, solidExp.Current() );
5341           }
5342         }
5343       }
5344       // make other sub-shapes computed
5345       setSubmeshesComputed( theMesh, theShape );
5346     }
5347
5348     // remove free nodes
5349     //if ( SMESHDS_SubMesh * smDS = meshDS->MeshElements( helper.GetSubShapeID() ))
5350     {
5351       std::vector< const SMDS_MeshNode* > nodesToRemove;
5352       // get intersection nodes
5353       for ( int iDir = 0; iDir < 3; ++iDir )
5354       {
5355         vector< GridLine >& lines = grid._lines[ iDir ];
5356         for ( size_t i = 0; i < lines.size(); ++i )
5357         {
5358           multiset< F_IntersectPoint >::iterator ip = lines[i]._intPoints.begin();
5359           for ( ; ip != lines[i]._intPoints.end(); ++ip )
5360             if ( ip->_node && ip->_node->NbInverseElements() == 0 && !ip->_node->isMarked() )
5361             {
5362               nodesToRemove.push_back( ip->_node );
5363               ip->_node->setIsMarked( true );
5364             }
5365         }
5366       }
5367       // get grid nodes
5368       for ( size_t i = 0; i < grid._nodes.size(); ++i )
5369         if ( grid._nodes[i] && grid._nodes[i]->NbInverseElements() == 0 &&
5370              !grid._nodes[i]->isMarked() )
5371         {
5372           nodesToRemove.push_back( grid._nodes[i] );
5373           grid._nodes[i]->setIsMarked( true );
5374         }
5375
5376       // do remove
5377       for ( size_t i = 0; i < nodesToRemove.size(); ++i )
5378         meshDS->RemoveFreeNode( nodesToRemove[i], /*smD=*/0, /*fromGroups=*/false );
5379     }
5380
5381     return nbAdded;
5382
5383   }
5384   // SMESH_ComputeError is not caught at SMESH_submesh level for an unknown reason
5385   catch ( SMESH_ComputeError& e)
5386   {
5387     return error( SMESH_ComputeErrorPtr( new SMESH_ComputeError( e )));
5388   }
5389   return false;
5390 }
5391
5392 //=============================================================================
5393 /*!
5394  *  Evaluate
5395  */
5396 //=============================================================================
5397
5398 bool StdMeshers_Cartesian_3D::Evaluate(SMESH_Mesh &         theMesh,
5399                                        const TopoDS_Shape & theShape,
5400                                        MapShapeNbElems&     theResMap)
5401 {
5402   // TODO
5403 //   std::vector<int> aResVec(SMDSEntity_Last);
5404 //   for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
5405 //   if(IsQuadratic) {
5406 //     aResVec[SMDSEntity_Quad_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
5407 //     int nb1d_face0_int = ( nb2d_face0*4 - nb1d ) / 2;
5408 //     aResVec[SMDSEntity_Node] = nb0d_face0 * ( 2*nb2d/nb1d - 1 ) - nb1d_face0_int * nb2d/nb1d;
5409 //   }
5410 //   else {
5411 //     aResVec[SMDSEntity_Node] = nb0d_face0 * ( nb2d/nb1d - 1 );
5412 //     aResVec[SMDSEntity_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
5413 //   }
5414 //   SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
5415 //   aResMap.insert(std::make_pair(sm,aResVec));
5416
5417   return true;
5418 }
5419
5420 //=============================================================================
5421 namespace
5422 {
5423   /*!
5424    * \brief Event listener setting/unsetting _alwaysComputed flag to
5425    *        submeshes of inferior levels to prevent their computing
5426    */
5427   struct _EventListener : public SMESH_subMeshEventListener
5428   {
5429     string _algoName;
5430
5431     _EventListener(const string& algoName):
5432       SMESH_subMeshEventListener(/*isDeletable=*/true,"StdMeshers_Cartesian_3D::_EventListener"),
5433       _algoName(algoName)
5434     {}
5435     // --------------------------------------------------------------------------------
5436     // setting/unsetting _alwaysComputed flag to submeshes of inferior levels
5437     //
5438     static void setAlwaysComputed( const bool     isComputed,
5439                                    SMESH_subMesh* subMeshOfSolid)
5440     {
5441       SMESH_subMeshIteratorPtr smIt =
5442         subMeshOfSolid->getDependsOnIterator(/*includeSelf=*/false, /*complexShapeFirst=*/false);
5443       while ( smIt->more() )
5444       {
5445         SMESH_subMesh* sm = smIt->next();
5446         sm->SetIsAlwaysComputed( isComputed );
5447       }
5448       subMeshOfSolid->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
5449     }
5450
5451     // --------------------------------------------------------------------------------
5452     // unsetting _alwaysComputed flag if "Cartesian_3D" was removed
5453     //
5454     virtual void ProcessEvent(const int          event,
5455                               const int          eventType,
5456                               SMESH_subMesh*     subMeshOfSolid,
5457                               SMESH_subMeshEventListenerData* data,
5458                               const SMESH_Hypothesis*         hyp = 0)
5459     {
5460       if ( eventType == SMESH_subMesh::COMPUTE_EVENT )
5461       {
5462         setAlwaysComputed( subMeshOfSolid->GetComputeState() == SMESH_subMesh::COMPUTE_OK,
5463                            subMeshOfSolid );
5464       }
5465       else
5466       {
5467         SMESH_Algo* algo3D = subMeshOfSolid->GetAlgo();
5468         if ( !algo3D || _algoName != algo3D->GetName() )
5469           setAlwaysComputed( false, subMeshOfSolid );
5470       }
5471     }
5472
5473     // --------------------------------------------------------------------------------
5474     // set the event listener
5475     //
5476     static void SetOn( SMESH_subMesh* subMeshOfSolid, const string& algoName )
5477     {
5478       subMeshOfSolid->SetEventListener( new _EventListener( algoName ),
5479                                         /*data=*/0,
5480                                         subMeshOfSolid );
5481     }
5482
5483   }; // struct _EventListener
5484
5485 } // namespace
5486
5487 //================================================================================
5488 /*!
5489  * \brief Sets event listener to submeshes if necessary
5490  *  \param subMesh - submesh where algo is set
5491  * This method is called when a submesh gets HYP_OK algo_state.
5492  * After being set, event listener is notified on each event of a submesh.
5493  */
5494 //================================================================================
5495
5496 void StdMeshers_Cartesian_3D::SetEventListener(SMESH_subMesh* subMesh)
5497 {
5498   _EventListener::SetOn( subMesh, GetName() );
5499 }
5500
5501 //================================================================================
5502 /*!
5503  * \brief Set _alwaysComputed flag to submeshes of inferior levels to avoid their computing
5504  */
5505 //================================================================================
5506
5507 void StdMeshers_Cartesian_3D::setSubmeshesComputed(SMESH_Mesh&         theMesh,
5508                                                    const TopoDS_Shape& theShape)
5509 {
5510   for ( TopExp_Explorer soExp( theShape, TopAbs_SOLID ); soExp.More(); soExp.Next() )
5511     _EventListener::setAlwaysComputed( true, theMesh.GetSubMesh( soExp.Current() ));
5512 }