Salome HOME
Fix regression of cartesian_algo.py example
[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     size_t getSolids( TGeomID ids[] );
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   size_t Hexahedron::getSolids( TGeomID ids[] )
2178   {
2179     if ( _grid->_geometry.IsOneSolid() )
2180     {
2181       ids[0] = _grid->GetSolid()->ID();
2182       return 1;
2183     }
2184     // count intersection points belonging to each SOLID
2185     TID2Nb id2NbPoints;
2186     id2NbPoints.reserve( 3 );
2187
2188     _origNodeInd = _grid->NodeIndex( _i,_j,_k );
2189     for ( int iN = 0; iN < 8; ++iN )
2190     {
2191       _hexNodes[iN]._node     = _grid->_nodes   [ _origNodeInd + _nodeShift[iN] ];
2192       _hexNodes[iN]._intPoint = _grid->_gridIntP[ _origNodeInd + _nodeShift[iN] ];
2193
2194       if ( _hexNodes[iN]._intPoint ) // intersection with a FACE
2195       {
2196         for ( size_t iF = 0; iF < _hexNodes[iN]._intPoint->_faceIDs.size(); ++iF )
2197         {
2198           const vector< TGeomID > & solidIDs =
2199             _grid->GetSolidIDs( _hexNodes[iN]._intPoint->_faceIDs[iF] );
2200           for ( size_t i = 0; i < solidIDs.size(); ++i )
2201             insertAndIncrement( solidIDs[i], id2NbPoints );
2202         }
2203       }
2204       else if ( _hexNodes[iN]._node ) // node inside a SOLID
2205       {
2206         insertAndIncrement( _hexNodes[iN]._node->GetShapeID(), id2NbPoints );
2207       }
2208     }
2209
2210     for ( int iL = 0; iL < 12; ++iL )
2211     {
2212       const _Link& link = _hexLinks[ iL ];
2213       for ( size_t iP = 0; iP < link._fIntPoints.size(); ++iP )
2214       {
2215         for ( size_t iF = 0; iF < link._fIntPoints[iP]->_faceIDs.size(); ++iF )
2216         {
2217           const vector< TGeomID > & solidIDs =
2218             _grid->GetSolidIDs( link._fIntPoints[iP]->_faceIDs[iF] );
2219           for ( size_t i = 0; i < solidIDs.size(); ++i )
2220             insertAndIncrement( solidIDs[i], id2NbPoints );
2221         }
2222       }
2223     }
2224
2225     for ( size_t iP = 0; iP < _eIntPoints.size(); ++iP )
2226     {
2227       const vector< TGeomID > & solidIDs = _grid->GetSolidIDs( _eIntPoints[iP]->_shapeID );
2228       for ( size_t i = 0; i < solidIDs.size(); ++i )
2229         insertAndIncrement( solidIDs[i], id2NbPoints );
2230     }
2231
2232     size_t nbSolids = 0;
2233     for ( TID2Nb::iterator id2nb = id2NbPoints.begin(); id2nb != id2NbPoints.end(); ++id2nb )
2234       if ( id2nb->second >= 3 )
2235         ids[ nbSolids++ ] = id2nb->first;
2236
2237     return nbSolids;
2238   }
2239
2240   //================================================================================
2241   /*!
2242    * \brief Count cuts by INTERNAL FACEs and set _Node::_isInternalFlags
2243    */
2244   bool Hexahedron::isCutByInternalFace( IsInternalFlag & maxFlag )
2245   {
2246     TID2Nb id2NbPoints;
2247     id2NbPoints.reserve( 3 );
2248
2249     for ( size_t iN = 0; iN < _intNodes.size(); ++iN )
2250       for ( size_t iF = 0; iF < _intNodes[iN]._intPoint->_faceIDs.size(); ++iF )
2251       {
2252         if ( _grid->IsInternal( _intNodes[iN]._intPoint->_faceIDs[iF]))
2253           insertAndIncrement( _intNodes[iN]._intPoint->_faceIDs[iF], id2NbPoints );
2254       }
2255     for ( size_t iN = 0; iN < 8; ++iN )
2256       if ( _hexNodes[iN]._intPoint )
2257         for ( size_t iF = 0; iF < _hexNodes[iN]._intPoint->_faceIDs.size(); ++iF )
2258         {
2259           if ( _grid->IsInternal( _hexNodes[iN]._intPoint->_faceIDs[iF]))
2260             insertAndIncrement( _hexNodes[iN]._intPoint->_faceIDs[iF], id2NbPoints );
2261         }
2262
2263     maxFlag = IS_NOT_INTERNAL;
2264     for ( TID2Nb::iterator id2nb = id2NbPoints.begin(); id2nb != id2NbPoints.end(); ++id2nb )
2265     {
2266       TGeomID        intFace = id2nb->first;
2267       IsInternalFlag intFlag = ( id2nb->second >= 3 ? IS_CUT_BY_INTERNAL_FACE : IS_INTERNAL );
2268       if ( intFlag > maxFlag )
2269         maxFlag = intFlag;
2270
2271       for ( size_t iN = 0; iN < _intNodes.size(); ++iN )
2272         if ( _intNodes[iN].IsOnFace( intFace ))
2273           _intNodes[iN].SetInternal( intFlag );
2274
2275       for ( size_t iN = 0; iN < 8; ++iN )
2276         if ( _hexNodes[iN].IsOnFace( intFace ))
2277           _hexNodes[iN].SetInternal( intFlag );
2278     }
2279
2280     return maxFlag;
2281   }
2282
2283   //================================================================================
2284   /*!
2285    * \brief Return any FACE interfering with this Hexahedron
2286    */
2287   TGeomID Hexahedron::getAnyFace() const
2288   {
2289     TID2Nb id2NbPoints;
2290     id2NbPoints.reserve( 3 );
2291
2292     for ( size_t iN = 0; iN < _intNodes.size(); ++iN )
2293       for ( size_t iF = 0; iF < _intNodes[iN]._intPoint->_faceIDs.size(); ++iF )
2294         insertAndIncrement( _intNodes[iN]._intPoint->_faceIDs[iF], id2NbPoints );
2295
2296     for ( size_t iN = 0; iN < 8; ++iN )
2297       if ( _hexNodes[iN]._intPoint )
2298         for ( size_t iF = 0; iF < _hexNodes[iN]._intPoint->_faceIDs.size(); ++iF )
2299           insertAndIncrement( _hexNodes[iN]._intPoint->_faceIDs[iF], id2NbPoints );
2300
2301     for ( unsigned int minNb = 3; minNb > 0; --minNb )
2302       for ( TID2Nb::iterator id2nb = id2NbPoints.begin(); id2nb != id2NbPoints.end(); ++id2nb )
2303         if ( id2nb->second >= minNb )
2304           return id2nb->first;
2305
2306     return 0;
2307   }
2308
2309   //================================================================================
2310   /*!
2311    * \brief Initializes IJK by Hexahedron index
2312    */
2313   void Hexahedron::setIJK( size_t iCell )
2314   {
2315     size_t iNbCell = _grid->_coords[0].size() - 1;
2316     size_t jNbCell = _grid->_coords[1].size() - 1;
2317     _i = iCell % iNbCell;
2318     _j = ( iCell % ( iNbCell * jNbCell )) / iNbCell;
2319     _k = iCell / iNbCell / jNbCell;
2320   }
2321
2322   //================================================================================
2323   /*!
2324    * \brief Initializes its data by given grid cell (countered from zero)
2325    */
2326   void Hexahedron::init( size_t iCell )
2327   {
2328     setIJK( iCell );
2329     init( _i, _j, _k );
2330   }
2331
2332   //================================================================================
2333   /*!
2334    * \brief Initializes its data by given grid cell nodes and intersections
2335    */
2336   void Hexahedron::init( size_t i, size_t j, size_t k, const Solid* solid )
2337   {
2338     _i = i; _j = j; _k = k;
2339
2340     if ( !solid )
2341       solid = _grid->GetSolid();
2342
2343     // set nodes of grid to nodes of the hexahedron and
2344     // count nodes at hexahedron corners located IN and ON geometry
2345     _nbCornerNodes = _nbBndNodes = 0;
2346     _origNodeInd   = _grid->NodeIndex( i,j,k );
2347     for ( int iN = 0; iN < 8; ++iN )
2348     {
2349       _hexNodes[iN]._isInternalFlags = 0;
2350
2351       _hexNodes[iN]._node     = _grid->_nodes   [ _origNodeInd + _nodeShift[iN] ];
2352       _hexNodes[iN]._intPoint = _grid->_gridIntP[ _origNodeInd + _nodeShift[iN] ];
2353
2354       if ( _hexNodes[iN]._node && !solid->Contains( _hexNodes[iN]._node->GetShapeID() ))
2355         _hexNodes[iN]._node = 0;
2356       if ( _hexNodes[iN]._intPoint && !solid->ContainsAny( _hexNodes[iN]._intPoint->_faceIDs ))
2357         _hexNodes[iN]._intPoint = 0;
2358
2359       _nbCornerNodes += bool( _hexNodes[iN]._node );
2360       _nbBndNodes    += bool( _hexNodes[iN]._intPoint );
2361     }
2362     _sideLength[0] = _grid->_coords[0][i+1] - _grid->_coords[0][i];
2363     _sideLength[1] = _grid->_coords[1][j+1] - _grid->_coords[1][j];
2364     _sideLength[2] = _grid->_coords[2][k+1] - _grid->_coords[2][k];
2365
2366     _intNodes.clear();
2367     _vIntNodes.clear();
2368
2369     if ( _nbFaceIntNodes + _eIntPoints.size()                  > 0 &&
2370          _nbFaceIntNodes + _eIntPoints.size() + _nbCornerNodes > 3)
2371     {
2372       _intNodes.reserve( 3 * _nbBndNodes + _nbFaceIntNodes + _eIntPoints.size() );
2373
2374       // this method can be called in parallel, so use own helper
2375       SMESH_MesherHelper helper( *_grid->_helper->GetMesh() );
2376
2377       // Create sub-links (_Link::_splits) by splitting links with _Link::_fIntPoints
2378       // ---------------------------------------------------------------
2379       _Link split;
2380       for ( int iLink = 0; iLink < 12; ++iLink )
2381       {
2382         _Link& link = _hexLinks[ iLink ];
2383         link._fIntNodes.clear();
2384         link._fIntNodes.reserve( link._fIntPoints.size() );
2385         for ( size_t i = 0; i < link._fIntPoints.size(); ++i )
2386           if ( solid->ContainsAny( link._fIntPoints[i]->_faceIDs ))
2387           {
2388             _intNodes.push_back( _Node( 0, link._fIntPoints[i] ));
2389             link._fIntNodes.push_back( & _intNodes.back() );
2390           }
2391
2392         link._splits.clear();
2393         split._nodes[ 0 ] = link._nodes[0];
2394         bool isOut = ( ! link._nodes[0]->Node() );
2395         bool checkTransition;
2396         for ( size_t i = 0; i < link._fIntNodes.size(); ++i )
2397         {
2398           const bool isGridNode = ( ! link._fIntNodes[i]->Node() );
2399           if ( !isGridNode ) // intersection non-coincident with a grid node
2400           {
2401             if ( split._nodes[ 0 ]->Node() && !isOut )
2402             {
2403               split._nodes[ 1 ] = link._fIntNodes[i];
2404               link._splits.push_back( split );
2405             }
2406             split._nodes[ 0 ] = link._fIntNodes[i];
2407             checkTransition = true;
2408           }
2409           else // FACE intersection coincident with a grid node (at link ends)
2410           {
2411             checkTransition = ( i == 0 && link._nodes[0]->Node() );
2412           }
2413           if ( checkTransition )
2414           {
2415             const vector< TGeomID >& faceIDs = link._fIntNodes[i]->_intPoint->_faceIDs;
2416             if ( _grid->IsInternal( faceIDs.back() ))
2417               isOut = false;
2418             else if ( faceIDs.size() > 1 || _eIntPoints.size() > 0 )
2419               isOut = isOutPoint( link, i, helper, solid );
2420             else
2421             {
2422               bool okTransi = _grid->IsCorrectTransition( faceIDs[0], solid );
2423               switch ( link._fIntNodes[i]->FaceIntPnt()->_transition ) {
2424               case Trans_OUT: isOut = okTransi;  break;
2425               case Trans_IN : isOut = !okTransi; break;
2426               default:
2427                 isOut = isOutPoint( link, i, helper, solid );
2428               }
2429             }
2430           }
2431         }
2432         if ( link._nodes[ 1 ]->Node() && split._nodes[ 0 ]->Node() && !isOut )
2433         {
2434           split._nodes[ 1 ] = link._nodes[1];
2435           link._splits.push_back( split );
2436         }
2437       }
2438
2439       // Create _Node's at intersections with EDGEs.
2440       // --------------------------------------------
2441       // 1) add this->_eIntPoints to _Face::_eIntNodes
2442       // 2) fill _intNodes and _vIntNodes
2443       //
2444       const double tol2 = _grid->_tol * _grid->_tol;
2445       int facets[3], nbFacets, subEntity;
2446
2447       for ( int iF = 0; iF < 6; ++iF )
2448         _hexQuads[ iF ]._eIntNodes.clear();
2449
2450       for ( size_t iP = 0; iP < _eIntPoints.size(); ++iP )
2451       {
2452         if ( !solid->ContainsAny( _eIntPoints[iP]->_faceIDs ))
2453           continue;
2454         nbFacets = getEntity( _eIntPoints[iP], facets, subEntity );
2455         _Node* equalNode = 0;
2456         switch( nbFacets ) {
2457         case 1: // in a _Face
2458         {
2459           _Face& quad = _hexQuads[ facets[0] - SMESH_Block::ID_FirstF ];
2460           equalNode = findEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
2461           if ( equalNode ) {
2462             equalNode->Add( _eIntPoints[ iP ] );
2463           }
2464           else {
2465             _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
2466             quad._eIntNodes.push_back( & _intNodes.back() );
2467           }
2468           break;
2469         }
2470         case 2: // on a _Link
2471         {
2472           _Link& link = _hexLinks[ subEntity - SMESH_Block::ID_FirstE ];
2473           if ( link._splits.size() > 0 )
2474           {
2475             equalNode = findEqualNode( link._fIntNodes, _eIntPoints[ iP ], tol2 );
2476             if ( equalNode )
2477               equalNode->Add( _eIntPoints[ iP ] );
2478             else if ( link._splits.size() == 1 &&
2479                       link._splits[0]._nodes[0] &&
2480                       link._splits[0]._nodes[1] )
2481               link._splits.clear(); // hex edge is divided by _eIntPoints[iP]
2482           }
2483           //else
2484           if ( !equalNode )
2485           {
2486             _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
2487             bool newNodeUsed = false;
2488             for ( int iF = 0; iF < 2; ++iF )
2489             {
2490               _Face& quad = _hexQuads[ facets[iF] - SMESH_Block::ID_FirstF ];
2491               equalNode = findEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
2492               if ( equalNode ) {
2493                 equalNode->Add( _eIntPoints[ iP ] );
2494               }
2495               else {
2496                 quad._eIntNodes.push_back( & _intNodes.back() );
2497                 newNodeUsed = true;
2498               }
2499             }
2500             if ( !newNodeUsed )
2501               _intNodes.pop_back();
2502           }
2503           break;
2504         }
2505         case 3: // at a corner
2506         {
2507           _Node& node = _hexNodes[ subEntity - SMESH_Block::ID_FirstV ];
2508           if ( node.Node() > 0 )
2509           {
2510             if ( node._intPoint )
2511               node._intPoint->Add( _eIntPoints[ iP ]->_faceIDs, _eIntPoints[ iP ]->_node );
2512           }
2513           else
2514           {
2515             _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
2516             for ( int iF = 0; iF < 3; ++iF )
2517             {
2518               _Face& quad = _hexQuads[ facets[iF] - SMESH_Block::ID_FirstF ];
2519               equalNode = findEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
2520               if ( equalNode ) {
2521                 equalNode->Add( _eIntPoints[ iP ] );
2522               }
2523               else {
2524                 quad._eIntNodes.push_back( & _intNodes.back() );
2525               }
2526             }
2527           }
2528           break;
2529         }
2530         } // switch( nbFacets )
2531
2532         if ( nbFacets == 0 ||
2533              _grid->ShapeType( _eIntPoints[ iP ]->_shapeID ) == TopAbs_VERTEX )
2534         {
2535           equalNode = findEqualNode( _vIntNodes, _eIntPoints[ iP ], tol2 );
2536           if ( equalNode ) {
2537             equalNode->Add( _eIntPoints[ iP ] );
2538           }
2539           else if ( nbFacets == 0 ) {
2540             if ( _intNodes.empty() || _intNodes.back().EdgeIntPnt() != _eIntPoints[ iP ])
2541               _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
2542             _vIntNodes.push_back( & _intNodes.back() );
2543           }
2544         }
2545       } // loop on _eIntPoints
2546     }
2547
2548     else if ( 3 < _nbCornerNodes && _nbCornerNodes < 8 ) // _nbFaceIntNodes == 0
2549     {
2550       _Link split;
2551       // create sub-links (_splits) of whole links
2552       for ( int iLink = 0; iLink < 12; ++iLink )
2553       {
2554         _Link& link = _hexLinks[ iLink ];
2555         link._splits.clear();
2556         if ( link._nodes[ 0 ]->Node() && link._nodes[ 1 ]->Node() )
2557         {
2558           split._nodes[ 0 ] = link._nodes[0];
2559           split._nodes[ 1 ] = link._nodes[1];
2560           link._splits.push_back( split );
2561         }
2562       }
2563     }
2564     return;
2565
2566   } // init( _i, _j, _k )
2567
2568   //================================================================================
2569   /*!
2570    * \brief Compute mesh volumes resulted from intersection of the Hexahedron
2571    */
2572   void Hexahedron::ComputeElements( const Solid* solid, int solidIndex )
2573   {
2574     if ( !solid )
2575     {
2576       solid = _grid->GetSolid();
2577       if ( !_grid->_geometry.IsOneSolid() )
2578       {
2579         TGeomID solidIDs[20];
2580         size_t nbSolids = getSolids( solidIDs );
2581         if ( nbSolids > 1 )
2582         {
2583           for ( size_t i = 0; i < nbSolids; ++i )
2584           {
2585             solid = _grid->GetSolid( solidIDs[i] );
2586             ComputeElements( solid, i );
2587             if ( !_volumeDefs._nodes.empty() && i < nbSolids - 1 )
2588               _volumeDefs.SetNext( new _volumeDef( _volumeDefs ));
2589           }
2590           return;
2591         }
2592         solid = _grid->GetSolid( solidIDs[0] );
2593       }
2594     }
2595
2596     init( _i, _j, _k, solid ); // get nodes and intersections from grid nodes and split links
2597
2598     int nbIntersections = _nbFaceIntNodes + _eIntPoints.size();
2599     if ( _nbCornerNodes + nbIntersections < 4 )
2600       return;
2601
2602     if ( _nbBndNodes == _nbCornerNodes && nbIntersections == 0 && isInHole() )
2603       return; // cell is in a hole
2604
2605     IsInternalFlag intFlag = IS_NOT_INTERNAL;
2606     if ( solid->HasInternalFaces() && this->isCutByInternalFace( intFlag ))
2607     {
2608       for ( _SplitIterator it( _hexLinks ); it.More(); it.Next() )
2609       {
2610         if ( compute( solid, intFlag ))
2611           _volumeDefs.SetNext( new _volumeDef( _volumeDefs ));
2612       }
2613     }
2614     else
2615     {
2616       if ( solidIndex >= 0 )
2617         intFlag = IS_CUT_BY_INTERNAL_FACE;
2618
2619       compute( solid, intFlag );
2620     }
2621   }
2622
2623   //================================================================================
2624   /*!
2625    * \brief Compute mesh volumes resulted from intersection of the Hexahedron
2626    */
2627   bool Hexahedron::compute( const Solid* solid, const IsInternalFlag intFlag )
2628   {
2629     _polygons.clear();
2630     _polygons.reserve( 20 );
2631
2632     for ( int iN = 0; iN < 8; ++iN )
2633       _hexNodes[iN]._usedInFace = 0;
2634
2635     // Create polygons from quadrangles
2636     // --------------------------------
2637
2638     vector< _OrientedLink > splits;
2639     vector<_Node*>          chainNodes;
2640     _Face*                  coplanarPolyg;
2641
2642     bool hasEdgeIntersections = !_eIntPoints.empty();
2643
2644     for ( int iF = 0; iF < 6; ++iF ) // loop on 6 sides of a hexahedron
2645     {
2646       _Face& quad = _hexQuads[ iF ] ;
2647
2648       _polygons.resize( _polygons.size() + 1 );
2649       _Face* polygon = &_polygons.back();
2650       polygon->_polyLinks.reserve( 20 );
2651
2652       splits.clear();
2653       for ( int iE = 0; iE < 4; ++iE ) // loop on 4 sides of a quadrangle
2654         for ( int iS = 0; iS < quad._links[ iE ].NbResultLinks(); ++iS )
2655           splits.push_back( quad._links[ iE ].ResultLink( iS ));
2656
2657       // add splits of links to a polygon and add _polyLinks to make
2658       // polygon's boundary closed
2659
2660       int nbSplits = splits.size();
2661       if (( nbSplits == 1 ) &&
2662           ( quad._eIntNodes.empty() ||
2663             splits[0].FirstNode()->IsLinked( splits[0].LastNode()->_intPoint )))
2664         //( quad._eIntNodes.empty() || _nbCornerNodes + nbIntersections > 6 ))
2665         nbSplits = 0;
2666
2667       for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
2668         if ( quad._eIntNodes[ iP ]->IsUsedInFace( polygon ))
2669           quad._eIntNodes[ iP ]->_usedInFace = 0;
2670
2671       size_t nbUsedEdgeNodes = 0;
2672       _Face* prevPolyg = 0; // polygon previously created from this quad
2673
2674       while ( nbSplits > 0 )
2675       {
2676         size_t iS = 0;
2677         while ( !splits[ iS ] )
2678           ++iS;
2679
2680         if ( !polygon->_links.empty() )
2681         {
2682           _polygons.resize( _polygons.size() + 1 );
2683           polygon = &_polygons.back();
2684           polygon->_polyLinks.reserve( 20 );
2685         }
2686         polygon->_links.push_back( splits[ iS ] );
2687         splits[ iS++ ]._link = 0;
2688         --nbSplits;
2689
2690         _Node* nFirst = polygon->_links.back().FirstNode();
2691         _Node *n1,*n2 = polygon->_links.back().LastNode();
2692         for ( ; nFirst != n2 && iS < splits.size(); ++iS )
2693         {
2694           _OrientedLink& split = splits[ iS ];
2695           if ( !split ) continue;
2696
2697           n1 = split.FirstNode();
2698           if ( n1 == n2 &&
2699                n1->_intPoint &&
2700                (( n1->_intPoint->_faceIDs.size() > 1 && isImplementEdges() ) ||
2701                 ( n1->_isInternalFlags )))
2702           {
2703             // n1 is at intersection with EDGE
2704             if ( findChainOnEdge( splits, polygon->_links.back(), split, iS, quad, chainNodes ))
2705             {
2706               for ( size_t i = 1; i < chainNodes.size(); ++i )
2707                 polygon->AddPolyLink( chainNodes[i-1], chainNodes[i], prevPolyg );
2708               if ( chainNodes.back() != n1 ) // not a partial cut by INTERNAL FACE
2709               {
2710                 prevPolyg = polygon;
2711                 n2 = chainNodes.back();
2712                 continue;
2713               }
2714             }
2715           }
2716           else if ( n1 != n2 )
2717           {
2718             // try to connect to intersections with EDGEs
2719             if ( quad._eIntNodes.size() > nbUsedEdgeNodes  &&
2720                  findChain( n2, n1, quad, chainNodes ))
2721             {
2722               for ( size_t i = 1; i < chainNodes.size(); ++i )
2723               {
2724                 polygon->AddPolyLink( chainNodes[i-1], chainNodes[i] );
2725                 nbUsedEdgeNodes += ( chainNodes[i]->IsUsedInFace( polygon ));
2726               }
2727               if ( chainNodes.back() != n1 )
2728               {
2729                 n2 = chainNodes.back();
2730                 --iS;
2731                 continue;
2732               }
2733             }
2734             // try to connect to a split ending on the same FACE
2735             else
2736             {
2737               _OrientedLink foundSplit;
2738               for ( size_t i = iS; i < splits.size() && !foundSplit; ++i )
2739                 if (( foundSplit = splits[ i ]) &&
2740                     ( n2->IsLinked( foundSplit.FirstNode()->_intPoint )))
2741                 {
2742                   iS = i - 1;
2743                 }
2744                 else
2745                 {
2746                   foundSplit._link = 0;
2747                 }
2748               if ( foundSplit )
2749               {
2750                 if ( n2 != foundSplit.FirstNode() )
2751                 {
2752                   polygon->AddPolyLink( n2, foundSplit.FirstNode() );
2753                   n2 = foundSplit.FirstNode();
2754                 }
2755                 continue;
2756               }
2757               else
2758               {
2759                 if ( n2->IsLinked( nFirst->_intPoint ))
2760                   break;
2761                 polygon->AddPolyLink( n2, n1, prevPolyg );
2762               }
2763             }
2764           } // if ( n1 != n2 )
2765
2766           polygon->_links.push_back( split );
2767           split._link = 0;
2768           --nbSplits;
2769           n2 = polygon->_links.back().LastNode();
2770
2771         } // loop on splits
2772
2773         if ( nFirst != n2 ) // close a polygon
2774         {
2775           if ( !findChain( n2, nFirst, quad, chainNodes ))
2776           {
2777             if ( !closePolygon( polygon, chainNodes ))
2778               if ( !isImplementEdges() )
2779                 chainNodes.push_back( nFirst );
2780           }
2781           for ( size_t i = 1; i < chainNodes.size(); ++i )
2782           {
2783             polygon->AddPolyLink( chainNodes[i-1], chainNodes[i], prevPolyg );
2784             nbUsedEdgeNodes += bool( chainNodes[i]->IsUsedInFace( polygon ));
2785           }
2786         }
2787
2788         if ( polygon->_links.size() < 3 && nbSplits > 0 )
2789         {
2790           polygon->_polyLinks.clear();
2791           polygon->_links.clear();
2792         }
2793       } // while ( nbSplits > 0 )
2794
2795       if ( polygon->_links.size() < 3 )
2796       {
2797         _polygons.pop_back();
2798       }
2799     }  // loop on 6 hexahedron sides
2800
2801     // Create polygons closing holes in a polyhedron
2802     // ----------------------------------------------
2803
2804     // clear _usedInFace
2805     for ( size_t iN = 0; iN < _intNodes.size(); ++iN )
2806       _intNodes[ iN ]._usedInFace = 0;
2807
2808     // add polygons to their links and mark used nodes
2809     for ( size_t iP = 0; iP < _polygons.size(); ++iP )
2810     {
2811       _Face& polygon = _polygons[ iP ];
2812       for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
2813       {
2814         polygon._links[ iL ].AddFace( &polygon );
2815         polygon._links[ iL ].FirstNode()->_usedInFace = &polygon;
2816       }
2817     }
2818     // find free links
2819     vector< _OrientedLink* > freeLinks;
2820     freeLinks.reserve(20);
2821     for ( size_t iP = 0; iP < _polygons.size(); ++iP )
2822     {
2823       _Face& polygon = _polygons[ iP ];
2824       for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
2825         if ( polygon._links[ iL ].NbFaces() < 2 )
2826           freeLinks.push_back( & polygon._links[ iL ]);
2827     }
2828     int nbFreeLinks = freeLinks.size();
2829     if ( nbFreeLinks == 1 ) return false;
2830
2831     // put not used intersection nodes to _vIntNodes
2832     int nbVertexNodes = 0; // nb not used vertex nodes
2833     {
2834       for ( size_t iN = 0; iN < _vIntNodes.size(); ++iN )
2835         nbVertexNodes += ( !_vIntNodes[ iN ]->IsUsedInFace() );
2836
2837       const double tol = 1e-3 * Min( Min( _sideLength[0], _sideLength[1] ), _sideLength[0] );
2838       for ( size_t iN = _nbFaceIntNodes; iN < _intNodes.size(); ++iN )
2839       {
2840         if ( _intNodes[ iN ].IsUsedInFace() ) continue;
2841         if ( dynamic_cast< const F_IntersectPoint* >( _intNodes[ iN ]._intPoint )) continue;
2842         _Node* equalNode =
2843           findEqualNode( _vIntNodes, _intNodes[ iN ].EdgeIntPnt(), tol*tol );
2844         if ( !equalNode )
2845         {
2846           _vIntNodes.push_back( &_intNodes[ iN ]);
2847           ++nbVertexNodes;
2848         }
2849       }
2850     }
2851
2852     set<TGeomID> usedFaceIDs;
2853     vector< TGeomID > faces;
2854     TGeomID curFace = 0;
2855     const size_t nbQuadPolygons = _polygons.size();
2856     E_IntersectPoint ipTmp;
2857
2858     // create polygons by making closed chains of free links
2859     size_t iPolygon = _polygons.size();
2860     while ( nbFreeLinks > 0 )
2861     {
2862       if ( iPolygon == _polygons.size() )
2863       {
2864         _polygons.resize( _polygons.size() + 1 );
2865         _polygons[ iPolygon ]._polyLinks.reserve( 20 );
2866         _polygons[ iPolygon ]._links.reserve( 20 );
2867       }
2868       _Face& polygon = _polygons[ iPolygon ];
2869
2870       _OrientedLink* curLink = 0;
2871       _Node*         curNode;
2872       if (( !hasEdgeIntersections ) ||
2873           ( nbFreeLinks < 4 && nbVertexNodes == 0 ))
2874       {
2875         // get a remaining link to start from
2876         for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2877           if (( curLink = freeLinks[ iL ] ))
2878             freeLinks[ iL ] = 0;
2879         polygon._links.push_back( *curLink );
2880         --nbFreeLinks;
2881         do
2882         {
2883           // find all links connected to curLink
2884           curNode = curLink->FirstNode();
2885           curLink = 0;
2886           for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2887             if ( freeLinks[ iL ] && freeLinks[ iL ]->LastNode() == curNode )
2888             {
2889               curLink = freeLinks[ iL ];
2890               freeLinks[ iL ] = 0;
2891               --nbFreeLinks;
2892               polygon._links.push_back( *curLink );
2893             }
2894         } while ( curLink );
2895       }
2896       else // there are intersections with EDGEs
2897       {
2898         // get a remaining link to start from, one lying on minimal nb of FACEs
2899         {
2900           typedef pair< TGeomID, int > TFaceOfLink;
2901           TFaceOfLink faceOfLink( -1, -1 );
2902           TFaceOfLink facesOfLink[3] = { faceOfLink, faceOfLink, faceOfLink };
2903           for ( size_t iL = 0; iL < freeLinks.size(); ++iL )
2904             if ( freeLinks[ iL ] )
2905             {
2906               faces = freeLinks[ iL ]->GetNotUsedFace( usedFaceIDs );
2907               if ( faces.size() == 1 )
2908               {
2909                 faceOfLink = TFaceOfLink( faces[0], iL );
2910                 if ( !freeLinks[ iL ]->HasEdgeNodes() )
2911                   break;
2912                 facesOfLink[0] = faceOfLink;
2913               }
2914               else if ( facesOfLink[0].first < 0 )
2915               {
2916                 faceOfLink = TFaceOfLink(( faces.empty() ? -1 : faces[0]), iL );
2917                 facesOfLink[ 1 + faces.empty() ] = faceOfLink;
2918               }
2919             }
2920           for ( int i = 0; faceOfLink.first < 0 && i < 3; ++i )
2921             faceOfLink = facesOfLink[i];
2922
2923           if ( faceOfLink.first < 0 ) // all faces used
2924           {
2925             for ( size_t iL = 0; iL < freeLinks.size() && faceOfLink.first < 1; ++iL )
2926               if (( curLink = freeLinks[ iL ]))
2927               {
2928                 faceOfLink.first = 
2929                   curLink->FirstNode()->IsLinked( curLink->LastNode()->_intPoint );
2930                 faceOfLink.second = iL;
2931               }
2932             usedFaceIDs.clear();
2933           }
2934           curFace = faceOfLink.first;
2935           curLink = freeLinks[ faceOfLink.second ];
2936           freeLinks[ faceOfLink.second ] = 0;
2937         }
2938         usedFaceIDs.insert( curFace );
2939         polygon._links.push_back( *curLink );
2940         --nbFreeLinks;
2941
2942         // find all links lying on a curFace
2943         do
2944         {
2945           // go forward from curLink
2946           curNode = curLink->LastNode();
2947           curLink = 0;
2948           for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2949             if ( freeLinks[ iL ] &&
2950                  freeLinks[ iL ]->FirstNode() == curNode &&
2951                  freeLinks[ iL ]->LastNode()->IsOnFace( curFace ))
2952             {
2953               curLink = freeLinks[ iL ];
2954               freeLinks[ iL ] = 0;
2955               polygon._links.push_back( *curLink );
2956               --nbFreeLinks;
2957             }
2958         } while ( curLink );
2959
2960         std::reverse( polygon._links.begin(), polygon._links.end() );
2961
2962         curLink = & polygon._links.back();
2963         do
2964         {
2965           // go backward from curLink
2966           curNode = curLink->FirstNode();
2967           curLink = 0;
2968           for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2969             if ( freeLinks[ iL ] &&
2970                  freeLinks[ iL ]->LastNode() == curNode &&
2971                  freeLinks[ iL ]->FirstNode()->IsOnFace( curFace ))
2972             {
2973               curLink = freeLinks[ iL ];
2974               freeLinks[ iL ] = 0;
2975               polygon._links.push_back( *curLink );
2976               --nbFreeLinks;
2977             }
2978         } while ( curLink );
2979
2980         curNode = polygon._links.back().FirstNode();
2981
2982         if ( polygon._links[0].LastNode() != curNode )
2983         {
2984           if ( nbVertexNodes > 0 )
2985           {
2986             // add links with _vIntNodes if not already used
2987             chainNodes.clear();
2988             for ( size_t iN = 0; iN < _vIntNodes.size(); ++iN )
2989               if ( !_vIntNodes[ iN ]->IsUsedInFace() &&
2990                    _vIntNodes[ iN ]->IsOnFace( curFace ))
2991               {
2992                 _vIntNodes[ iN ]->_usedInFace = &polygon;
2993                 chainNodes.push_back( _vIntNodes[ iN ] );
2994               }
2995             if ( chainNodes.size() > 1 &&
2996                  curFace != _grid->PseudoIntExtFaceID() ) /////// TODO
2997             {
2998               sortVertexNodes( chainNodes, curNode, curFace );
2999             }
3000             for ( size_t i = 0; i < chainNodes.size(); ++i )
3001             {
3002               polygon.AddPolyLink( chainNodes[ i ], curNode );
3003               curNode = chainNodes[ i ];
3004               freeLinks.push_back( &polygon._links.back() );
3005               ++nbFreeLinks;
3006             }
3007             nbVertexNodes -= chainNodes.size();
3008           }
3009           // if ( polygon._links.size() > 1 )
3010           {
3011             polygon.AddPolyLink( polygon._links[0].LastNode(), curNode );
3012             freeLinks.push_back( &polygon._links.back() );
3013             ++nbFreeLinks;
3014           }
3015         }
3016       } // if there are intersections with EDGEs
3017
3018       if ( polygon._links.size() < 2 ||
3019            polygon._links[0].LastNode() != polygon._links.back().FirstNode() )
3020         return false; // closed polygon not found -> invalid polyhedron
3021
3022       if ( polygon._links.size() == 2 )
3023       {
3024         if ( freeLinks.back() == &polygon._links.back() )
3025         {
3026           freeLinks.pop_back();
3027           --nbFreeLinks;
3028         }
3029         if ( polygon._links.front().NbFaces() > 0 )
3030           polygon._links.back().AddFace( polygon._links.front()._link->_faces[0] );
3031         if ( polygon._links.back().NbFaces() > 0 )
3032           polygon._links.front().AddFace( polygon._links.back()._link->_faces[0] );
3033
3034         if ( iPolygon == _polygons.size()-1 )
3035           _polygons.pop_back();
3036       }
3037       else // polygon._links.size() >= 2
3038       {
3039         // add polygon to its links
3040         for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
3041         {
3042           polygon._links[ iL ].AddFace( &polygon );
3043           polygon._links[ iL ].Reverse();
3044         }
3045         if ( /*hasEdgeIntersections &&*/ iPolygon == _polygons.size() - 1 )
3046         {
3047           // check that a polygon does not lie on a hexa side
3048           coplanarPolyg = 0;
3049           for ( size_t iL = 0; iL < polygon._links.size() && !coplanarPolyg; ++iL )
3050           {
3051             if ( polygon._links[ iL ].NbFaces() < 2 )
3052               continue; // it's a just added free link
3053             // look for a polygon made on a hexa side and sharing
3054             // two or more haxa links
3055             size_t iL2;
3056             coplanarPolyg = polygon._links[ iL ]._link->_faces[0];
3057             for ( iL2 = iL + 1; iL2 < polygon._links.size(); ++iL2 )
3058               if ( polygon._links[ iL2 ]._link->_faces[0] == coplanarPolyg &&
3059                    !coplanarPolyg->IsPolyLink( polygon._links[ iL  ]) &&
3060                    !coplanarPolyg->IsPolyLink( polygon._links[ iL2 ]) &&
3061                    coplanarPolyg < & _polygons[ nbQuadPolygons ])
3062                 break;
3063             if ( iL2 == polygon._links.size() )
3064               coplanarPolyg = 0;
3065           }
3066           if ( coplanarPolyg ) // coplanar polygon found
3067           {
3068             freeLinks.resize( freeLinks.size() - polygon._polyLinks.size() );
3069             nbFreeLinks -= polygon._polyLinks.size();
3070
3071             // an E_IntersectPoint used to mark nodes of coplanarPolyg
3072             // as lying on curFace while they are not at intersection with geometry
3073             ipTmp._faceIDs.resize(1);
3074             ipTmp._faceIDs[0] = curFace;
3075
3076             // fill freeLinks with links not shared by coplanarPolyg and polygon
3077             for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
3078               if ( polygon._links[ iL ]._link->_faces[1] &&
3079                    polygon._links[ iL ]._link->_faces[0] != coplanarPolyg )
3080               {
3081                 _Face* p = polygon._links[ iL ]._link->_faces[0];
3082                 for ( size_t iL2 = 0; iL2 < p->_links.size(); ++iL2 )
3083                   if ( p->_links[ iL2 ]._link == polygon._links[ iL ]._link )
3084                   {
3085                     freeLinks.push_back( & p->_links[ iL2 ] );
3086                     ++nbFreeLinks;
3087                     freeLinks.back()->RemoveFace( &polygon );
3088                     break;
3089                   }
3090               }
3091             for ( size_t iL = 0; iL < coplanarPolyg->_links.size(); ++iL )
3092               if ( coplanarPolyg->_links[ iL ]._link->_faces[1] &&
3093                    coplanarPolyg->_links[ iL ]._link->_faces[1] != &polygon )
3094               {
3095                 _Face* p = coplanarPolyg->_links[ iL ]._link->_faces[0];
3096                 if ( p == coplanarPolyg )
3097                   p = coplanarPolyg->_links[ iL ]._link->_faces[1];
3098                 for ( size_t iL2 = 0; iL2 < p->_links.size(); ++iL2 )
3099                   if ( p->_links[ iL2 ]._link == coplanarPolyg->_links[ iL ]._link )
3100                   {
3101                     // set links of coplanarPolyg in place of used freeLinks
3102                     // to re-create coplanarPolyg next
3103                     size_t iL3 = 0;
3104                     for ( ; iL3 < freeLinks.size() && freeLinks[ iL3 ]; ++iL3 );
3105                     if ( iL3 < freeLinks.size() )
3106                       freeLinks[ iL3 ] = ( & p->_links[ iL2 ] );
3107                     else
3108                       freeLinks.push_back( & p->_links[ iL2 ] );
3109                     ++nbFreeLinks;
3110                     freeLinks[ iL3 ]->RemoveFace( coplanarPolyg );
3111                     //  mark nodes of coplanarPolyg as lying on curFace
3112                     for ( int iN = 0; iN < 2; ++iN )
3113                     {
3114                       _Node* n = freeLinks[ iL3 ]->_link->_nodes[ iN ];
3115                       if ( n->_intPoint ) n->_intPoint->Add( ipTmp._faceIDs );
3116                       else                n->_intPoint = &ipTmp;
3117                     }
3118                     break;
3119                   }
3120               }
3121             // set coplanarPolyg to be re-created next
3122             for ( size_t iP = 0; iP < _polygons.size(); ++iP )
3123               if ( coplanarPolyg == & _polygons[ iP ] )
3124               {
3125                 iPolygon = iP;
3126                 _polygons[ iPolygon ]._links.clear();
3127                 _polygons[ iPolygon ]._polyLinks.clear();
3128                 break;
3129               }
3130             _polygons.pop_back();
3131             usedFaceIDs.erase( curFace );
3132             continue;
3133           } // if ( coplanarPolyg )
3134         } // if ( hasEdgeIntersections ) - search for coplanarPolyg
3135
3136         iPolygon = _polygons.size();
3137
3138       } // end of case ( polygon._links.size() > 2 )
3139     } // while ( nbFreeLinks > 0 )
3140
3141     // check volume size
3142     _hasTooSmall = ! checkPolyhedronSize( intFlag & IS_CUT_BY_INTERNAL_FACE );
3143
3144     for ( size_t i = 0; i < 8; ++i )
3145       if ( _hexNodes[ i ]._intPoint == &ipTmp )
3146         _hexNodes[ i ]._intPoint = 0;
3147
3148     if ( _hasTooSmall )
3149       return false; // too small volume
3150
3151     // create a classic cell if possible
3152
3153     int nbPolygons = 0;
3154     for ( size_t iF = 0; iF < _polygons.size(); ++iF )
3155       nbPolygons += (_polygons[ iF ]._links.size() > 0 );
3156
3157     //const int nbNodes = _nbCornerNodes + nbIntersections;
3158     int nbNodes = 0;
3159     for ( size_t i = 0; i < 8; ++i )
3160       nbNodes += _hexNodes[ i ].IsUsedInFace();
3161     for ( size_t i = 0; i < _intNodes.size(); ++i )
3162       nbNodes += _intNodes[ i ].IsUsedInFace();
3163
3164     bool isClassicElem = false;
3165     if (      nbNodes == 8 && nbPolygons == 6 ) isClassicElem = addHexa();
3166     else if ( nbNodes == 4 && nbPolygons == 4 ) isClassicElem = addTetra();
3167     else if ( nbNodes == 6 && nbPolygons == 5 ) isClassicElem = addPenta();
3168     else if ( nbNodes == 5 && nbPolygons == 5 ) isClassicElem = addPyra ();
3169     if ( !isClassicElem )
3170     {
3171       _volumeDefs._nodes.clear();
3172       _volumeDefs._quantities.clear();
3173
3174       for ( size_t iF = 0; iF < _polygons.size(); ++iF )
3175       {
3176         const size_t nbLinks = _polygons[ iF ]._links.size();
3177         if ( nbLinks == 0 ) continue;
3178         _volumeDefs._quantities.push_back( nbLinks );
3179         for ( size_t iL = 0; iL < nbLinks; ++iL )
3180           _volumeDefs._nodes.push_back( _polygons[ iF ]._links[ iL ].FirstNode() );
3181       }
3182     }
3183     _volumeDefs._solidID = solid->ID();
3184
3185     return !_volumeDefs._nodes.empty();
3186   }
3187   //================================================================================
3188   /*!
3189    * \brief Create elements in the mesh
3190    */
3191   int Hexahedron::MakeElements(SMESH_MesherHelper&                      helper,
3192                                const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap)
3193   {
3194     SMESHDS_Mesh* mesh = helper.GetMeshDS();
3195
3196     CellsAroundLink c( _grid, 0 );
3197     const size_t nbGridCells = c._nbCells[0] * c._nbCells[1] * c._nbCells[2];
3198     vector< Hexahedron* > allHexa( nbGridCells, 0 );
3199     int nbIntHex = 0;
3200
3201     // set intersection nodes from GridLine's to links of allHexa
3202     int i,j,k, cellIndex;
3203     for ( int iDir = 0; iDir < 3; ++iDir )
3204     {
3205       // loop on GridLine's parallel to iDir
3206       LineIndexer lineInd = _grid->GetLineIndexer( iDir );
3207       CellsAroundLink fourCells( _grid, iDir );
3208       for ( ; lineInd.More(); ++lineInd )
3209       {
3210         GridLine& line = _grid->_lines[ iDir ][ lineInd.LineIndex() ];
3211         multiset< F_IntersectPoint >::const_iterator ip = line._intPoints.begin();
3212         for ( ; ip != line._intPoints.end(); ++ip )
3213         {
3214           // if ( !ip->_node ) continue; // intersection at a grid node
3215           lineInd.SetIndexOnLine( ip->_indexOnLine );
3216           fourCells.Init( lineInd.I(), lineInd.J(), lineInd.K() );
3217           for ( int iL = 0; iL < 4; ++iL ) // loop on 4 cells sharing a link
3218           {
3219             if ( !fourCells.GetCell( iL, i,j,k, cellIndex ))
3220               continue;
3221             Hexahedron *& hex = allHexa[ cellIndex ];
3222             if ( !hex)
3223             {
3224               hex = new Hexahedron( *this, i, j, k, cellIndex );
3225               ++nbIntHex;
3226             }
3227             const int iLink = iL + iDir * 4;
3228             hex->_hexLinks[iLink]._fIntPoints.push_back( &(*ip) );
3229             hex->_nbFaceIntNodes += bool( ip->_node );
3230           }
3231         }
3232       }
3233     }
3234
3235     // implement geom edges into the mesh
3236     addEdges( helper, allHexa, edge2faceIDsMap );
3237
3238     // add not split hexahedra to the mesh
3239     int nbAdded = 0;
3240     TGeomID solidIDs[20];
3241     vector< Hexahedron* > intHexa; intHexa.reserve( nbIntHex );
3242     vector< const SMDS_MeshElement* > boundaryVolumes; boundaryVolumes.reserve( nbIntHex * 1.1 );
3243     for ( size_t i = 0; i < allHexa.size(); ++i )
3244     {
3245       // initialize this by not cut allHexa[ i ]
3246       Hexahedron * & hex = allHexa[ i ];
3247       if ( hex ) // split hexahedron
3248       {
3249         intHexa.push_back( hex );
3250         if ( hex->_nbFaceIntNodes > 0 || hex->_eIntPoints.size() > 0 )
3251           continue; // treat intersected hex later in parallel
3252         this->init( hex->_i, hex->_j, hex->_k );
3253       }
3254       else
3255       {
3256         this->init( i ); // == init(i,j,k)
3257       }
3258       if (( _nbCornerNodes == 8 ) &&
3259           ( _nbBndNodes < _nbCornerNodes || !isInHole() ))
3260       {
3261         // order of _hexNodes is defined by enum SMESH_Block::TShapeID
3262         SMDS_MeshElement* el =
3263           mesh->AddVolume( _hexNodes[0].Node(), _hexNodes[2].Node(),
3264                            _hexNodes[3].Node(), _hexNodes[1].Node(),
3265                            _hexNodes[4].Node(), _hexNodes[6].Node(),
3266                            _hexNodes[7].Node(), _hexNodes[5].Node() );
3267         TGeomID solidID = 0;
3268         if ( _nbBndNodes < _nbCornerNodes )
3269         {
3270           for ( int iN = 0; iN < 8 &&  !solidID; ++iN )
3271             if ( !_hexNodes[iN]._intPoint ) // no intersection
3272               solidID = _hexNodes[iN].Node()->GetShapeID();
3273         }
3274         else
3275         {
3276           getSolids( solidIDs );
3277           solidID = solidIDs[0];
3278         }
3279         mesh->SetMeshElementOnShape( el, solidID );
3280         ++nbAdded;
3281         if ( hex )
3282           intHexa.pop_back();
3283         if ( _grid->_toCreateFaces && _nbBndNodes >= 3 )
3284         {
3285           boundaryVolumes.push_back( el );
3286           el->setIsMarked( true );
3287         }
3288       }
3289       else if ( _nbCornerNodes > 3 && !hex )
3290       {
3291         // all intersection of hex with geometry are at grid nodes
3292         hex = new Hexahedron( *this, _i, _j, _k, i );
3293         intHexa.push_back( hex );
3294       }
3295     }
3296
3297     // add elements resulted from hexadron intersection
3298 #ifdef WITH_TBB
3299     tbb::parallel_for ( tbb::blocked_range<size_t>( 0, intHexa.size() ),
3300                         ParallelHexahedron( intHexa ),
3301                         tbb::simple_partitioner()); // ComputeElements() is called here
3302     for ( size_t i = 0; i < intHexa.size(); ++i )
3303       if ( Hexahedron * hex = intHexa[ i ] )
3304         nbAdded += hex->addVolumes( helper );
3305 #else
3306     for ( size_t i = 0; i < intHexa.size(); ++i )
3307       if ( Hexahedron * hex = intHexa[ i ] )
3308       {
3309         hex->ComputeElements();
3310         nbAdded += hex->addVolumes( helper );
3311       }
3312 #endif
3313
3314     // fill boundaryVolumes with volumes neighboring too small skipped volumes
3315     if ( _grid->_toCreateFaces )
3316     {
3317       for ( size_t i = 0; i < intHexa.size(); ++i )
3318         if ( Hexahedron * hex = intHexa[ i ] )
3319           hex->getBoundaryElems( boundaryVolumes );
3320     }
3321
3322     // create boundary mesh faces
3323     addFaces( helper, boundaryVolumes );
3324
3325     // create mesh edges
3326     addSegments( helper, edge2faceIDsMap );
3327
3328     for ( size_t i = 0; i < allHexa.size(); ++i )
3329       if ( allHexa[ i ] )
3330         delete allHexa[ i ];
3331
3332     return nbAdded;
3333   }
3334
3335   //================================================================================
3336   /*!
3337    * \brief Implements geom edges into the mesh
3338    */
3339   void Hexahedron::addEdges(SMESH_MesherHelper&                      helper,
3340                             vector< Hexahedron* >&                   hexes,
3341                             const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap)
3342   {
3343     if ( edge2faceIDsMap.empty() ) return;
3344
3345     // Prepare planes for intersecting with EDGEs
3346     GridPlanes pln[3];
3347     {
3348       for ( int iDirZ = 0; iDirZ < 3; ++iDirZ ) // iDirZ gives normal direction to planes
3349       {
3350         GridPlanes& planes = pln[ iDirZ ];
3351         int iDirX = ( iDirZ + 1 ) % 3;
3352         int iDirY = ( iDirZ + 2 ) % 3;
3353         planes._zNorm  = ( _grid->_axes[ iDirX ] ^ _grid->_axes[ iDirY ] ).Normalized();
3354         planes._zProjs.resize ( _grid->_coords[ iDirZ ].size() );
3355         planes._zProjs [0] = 0;
3356         const double       zFactor = _grid->_axes[ iDirZ ] * planes._zNorm;
3357         const vector< double > & u = _grid->_coords[ iDirZ ];
3358         for ( size_t i = 1; i < planes._zProjs.size(); ++i )
3359         {
3360           planes._zProjs [i] = zFactor * ( u[i] - u[0] );
3361         }
3362       }
3363     }
3364     const double deflection = _grid->_minCellSize / 20.;
3365     const double tol        = _grid->_tol;
3366     E_IntersectPoint ip;
3367
3368     TColStd_MapOfInteger intEdgeIDs; // IDs of not shared INTERNAL EDGES
3369
3370     // Intersect EDGEs with the planes
3371     map< TGeomID, vector< TGeomID > >::const_iterator e2fIt = edge2faceIDsMap.begin();
3372     for ( ; e2fIt != edge2faceIDsMap.end(); ++e2fIt )
3373     {
3374       const TGeomID  edgeID = e2fIt->first;
3375       const TopoDS_Edge & E = TopoDS::Edge( _grid->Shape( edgeID ));
3376       BRepAdaptor_Curve curve( E );
3377       TopoDS_Vertex v1 = helper.IthVertex( 0, E, false );
3378       TopoDS_Vertex v2 = helper.IthVertex( 1, E, false );
3379
3380       ip._faceIDs = e2fIt->second;
3381       ip._shapeID = edgeID;
3382
3383       bool isInternal = ( ip._faceIDs.size() == 1 && _grid->IsInternal( edgeID ));
3384       if ( isInternal )
3385       {
3386         intEdgeIDs.Add( edgeID );
3387         intEdgeIDs.Add( _grid->ShapeID( v1 ));
3388         intEdgeIDs.Add( _grid->ShapeID( v2 ));
3389       }
3390
3391       // discretize the EDGE
3392       GCPnts_UniformDeflection discret( curve, deflection, true );
3393       if ( !discret.IsDone() || discret.NbPoints() < 2 )
3394         continue;
3395
3396       // perform intersection
3397       E_IntersectPoint* eip, *vip;
3398       for ( int iDirZ = 0; iDirZ < 3; ++iDirZ )
3399       {
3400         GridPlanes& planes = pln[ iDirZ ];
3401         int      iDirX = ( iDirZ + 1 ) % 3;
3402         int      iDirY = ( iDirZ + 2 ) % 3;
3403         double    xLen = _grid->_coords[ iDirX ].back() - _grid->_coords[ iDirX ][0];
3404         double    yLen = _grid->_coords[ iDirY ].back() - _grid->_coords[ iDirY ][0];
3405         double    zLen = _grid->_coords[ iDirZ ].back() - _grid->_coords[ iDirZ ][0];
3406         int dIJK[3], d000[3] = { 0,0,0 };
3407         double o[3] = { _grid->_coords[0][0],
3408                         _grid->_coords[1][0],
3409                         _grid->_coords[2][0] };
3410
3411         // locate the 1st point of a segment within the grid
3412         gp_XYZ p1     = discret.Value( 1 ).XYZ();
3413         double u1     = discret.Parameter( 1 );
3414         double zProj1 = planes._zNorm * ( p1 - _grid->_origin );
3415
3416         _grid->ComputeUVW( p1, ip._uvw );
3417         int iX1 = int(( ip._uvw[iDirX] - o[iDirX]) / xLen * (_grid->_coords[ iDirX ].size() - 1));
3418         int iY1 = int(( ip._uvw[iDirY] - o[iDirY]) / yLen * (_grid->_coords[ iDirY ].size() - 1));
3419         int iZ1 = int(( ip._uvw[iDirZ] - o[iDirZ]) / zLen * (_grid->_coords[ iDirZ ].size() - 1));
3420         locateValue( iX1, ip._uvw[iDirX], _grid->_coords[ iDirX ], dIJK[ iDirX ], tol );
3421         locateValue( iY1, ip._uvw[iDirY], _grid->_coords[ iDirY ], dIJK[ iDirY ], tol );
3422         locateValue( iZ1, ip._uvw[iDirZ], _grid->_coords[ iDirZ ], dIJK[ iDirZ ], tol );
3423
3424         int ijk[3]; // grid index where a segment intersects a plane
3425         ijk[ iDirX ] = iX1;
3426         ijk[ iDirY ] = iY1;
3427         ijk[ iDirZ ] = iZ1;
3428
3429         // add the 1st vertex point to a hexahedron
3430         if ( iDirZ == 0 )
3431         {
3432           ip._point   = p1;
3433           ip._shapeID = _grid->ShapeID( v1 );
3434           vip = _grid->Add( ip );
3435           if ( isInternal )
3436             vip->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
3437           if ( !addIntersection( vip, hexes, ijk, d000 ))
3438             _grid->Remove( vip );
3439           ip._shapeID = edgeID;
3440         }
3441         for ( int iP = 2; iP <= discret.NbPoints(); ++iP )
3442         {
3443           // locate the 2nd point of a segment within the grid
3444           gp_XYZ p2     = discret.Value( iP ).XYZ();
3445           double u2     = discret.Parameter( iP );
3446           double zProj2 = planes._zNorm * ( p2 - _grid->_origin );
3447           int    iZ2    = iZ1;
3448           if ( Abs( zProj2 - zProj1 ) > std::numeric_limits<double>::min() )
3449           {
3450             locateValue( iZ2, zProj2, planes._zProjs, dIJK[ iDirZ ], tol );
3451
3452             // treat intersections with planes between 2 end points of a segment
3453             int dZ = ( iZ1 <= iZ2 ) ? +1 : -1;
3454             int iZ = iZ1 + ( iZ1 < iZ2 );
3455             for ( int i = 0, nb = Abs( iZ1 - iZ2 ); i < nb; ++i, iZ += dZ )
3456             {
3457               ip._point = findIntPoint( u1, zProj1, u2, zProj2,
3458                                         planes._zProjs[ iZ ],
3459                                         curve, planes._zNorm, _grid->_origin );
3460               _grid->ComputeUVW( ip._point.XYZ(), ip._uvw );
3461               locateValue( ijk[iDirX], ip._uvw[iDirX], _grid->_coords[iDirX], dIJK[iDirX], tol );
3462               locateValue( ijk[iDirY], ip._uvw[iDirY], _grid->_coords[iDirY], dIJK[iDirY], tol );
3463               ijk[ iDirZ ] = iZ;
3464
3465               // add ip to hex "above" the plane
3466               eip = _grid->Add( ip );
3467               if ( isInternal )
3468                 eip->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
3469               dIJK[ iDirZ ] = 0;
3470               bool added = addIntersection( eip, hexes, ijk, dIJK);
3471
3472               // add ip to hex "below" the plane
3473               ijk[ iDirZ ] = iZ-1;
3474               if ( !addIntersection( eip, hexes, ijk, dIJK ) &&
3475                    !added )
3476                 _grid->Remove( eip );
3477             }
3478           }
3479           iZ1    = iZ2;
3480           p1     = p2;
3481           u1     = u2;
3482           zProj1 = zProj2;
3483         }
3484         // add the 2nd vertex point to a hexahedron
3485         if ( iDirZ == 0 )
3486         {
3487           ip._point   = p1;
3488           ip._shapeID = _grid->ShapeID( v2 );
3489           _grid->ComputeUVW( p1, ip._uvw );
3490           locateValue( ijk[iDirX], ip._uvw[iDirX], _grid->_coords[iDirX], dIJK[iDirX], tol );
3491           locateValue( ijk[iDirY], ip._uvw[iDirY], _grid->_coords[iDirY], dIJK[iDirY], tol );
3492           ijk[ iDirZ ] = iZ1;
3493           bool sameV = ( v1.IsSame( v2 ));
3494           if ( !sameV )
3495             vip = _grid->Add( ip );
3496           if ( isInternal && !sameV )
3497             vip->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
3498           if ( !addIntersection( vip, hexes, ijk, d000 ) && !sameV )
3499             _grid->Remove( vip );
3500           ip._shapeID = edgeID;
3501         }
3502       } // loop on 3 grid directions
3503     } // loop on EDGEs
3504
3505
3506     if ( intEdgeIDs.Size() > 0 )
3507       cutByExtendedInternal( hexes, intEdgeIDs );
3508
3509     return;
3510   }
3511
3512   //================================================================================
3513   /*!
3514    * \brief Fully cut hexes that are partially cut by INTERNAL FACE.
3515    *        Cut them by extended INTERNAL FACE.
3516    */
3517   void Hexahedron::cutByExtendedInternal( std::vector< Hexahedron* >& hexes,
3518                                           const TColStd_MapOfInteger& intEdgeIDs )
3519   {
3520     IntAna_IntConicQuad intersection;
3521     SMESHDS_Mesh* meshDS = _grid->_helper->GetMeshDS();
3522     const double tol2 = _grid->_tol * _grid->_tol;
3523
3524     for ( size_t iH = 0; iH < hexes.size(); ++iH )
3525     {
3526       Hexahedron* hex = hexes[ iH ];
3527       if ( !hex || hex->_eIntPoints.size() < 2 )
3528         continue;
3529       if ( !intEdgeIDs.Contains( hex->_eIntPoints.back()->_shapeID ))
3530         continue;
3531
3532       // get 3 points on INTERNAL FACE to construct a cutting plane
3533       gp_Pnt p1 = hex->_eIntPoints[0]->_point;
3534       gp_Pnt p2 = hex->_eIntPoints[1]->_point;
3535       gp_Pnt p3 = hex->mostDistantInternalPnt( iH, p1, p2 );
3536
3537       gp_Vec norm = gp_Vec( p1, p2 ) ^ gp_Vec( p1, p3 );
3538       gp_Pln pln;
3539       try {
3540         pln = gp_Pln( p1, norm );
3541       }
3542       catch(...)
3543       {
3544         continue;
3545       }
3546
3547       TGeomID intFaceID = hex->_eIntPoints.back()->_faceIDs.front(); // FACE being "extended"
3548       TGeomID   solidID = _grid->GetSolid( intFaceID )->ID();
3549
3550       // cut links by the plane
3551       //bool isCut = false;
3552       for ( int iLink = 0; iLink < 12; ++iLink )
3553       {
3554         _Link& link = hex->_hexLinks[ iLink ];
3555         if ( !link._fIntPoints.empty() )
3556         {
3557           // if ( link._fIntPoints[0]->_faceIDs.back() == _grid->PseudoIntExtFaceID() )
3558           //   isCut = true;
3559           continue; // already cut link
3560         }
3561         if ( !link._nodes[0]->Node() ||
3562              !link._nodes[1]->Node() )
3563           continue; // outside link
3564
3565         if ( link._nodes[0]->IsOnFace( intFaceID ))
3566         {
3567           if ( link._nodes[0]->_intPoint->_faceIDs.back() != _grid->PseudoIntExtFaceID() )
3568             if ( p1.SquareDistance( link._nodes[0]->Point() ) < tol2  ||
3569                  p2.SquareDistance( link._nodes[0]->Point() ) < tol2 )
3570               link._nodes[0]->_intPoint->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
3571           continue; // link is cut by FACE being "extended"
3572         }
3573         if ( link._nodes[1]->IsOnFace( intFaceID ))
3574         {
3575           if ( link._nodes[1]->_intPoint->_faceIDs.back() != _grid->PseudoIntExtFaceID() )
3576             if ( p1.SquareDistance( link._nodes[1]->Point() ) < tol2  ||
3577                  p2.SquareDistance( link._nodes[1]->Point() ) < tol2 )
3578               link._nodes[1]->_intPoint->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
3579           continue; // link is cut by FACE being "extended"
3580         }
3581         gp_Pnt p4 = link._nodes[0]->Point();
3582         gp_Pnt p5 = link._nodes[1]->Point();
3583         gp_Lin line( p4, gp_Vec( p4, p5 ));
3584
3585         intersection.Perform( line, pln );
3586         if ( !intersection.IsDone() ||
3587              intersection.IsInQuadric() ||
3588              intersection.IsParallel() ||
3589              intersection.NbPoints() < 1 )
3590           continue;
3591
3592         double u = intersection.ParamOnConic(1);
3593         if ( u + _grid->_tol < 0 )
3594           continue;
3595         int       iDir = iLink / 4;
3596         int      index = (&hex->_i)[iDir];
3597         double linkLen = _grid->_coords[iDir][index+1] - _grid->_coords[iDir][index];
3598         if ( u - _grid->_tol > linkLen )
3599           continue;
3600
3601         if ( u < _grid->_tol ||
3602              u > linkLen - _grid->_tol ) // intersection at grid node
3603         {
3604           int  i = ! ( u < _grid->_tol ); // [0,1]
3605           int iN = link._nodes[ i ] - hex->_hexNodes; // [0-7]
3606
3607           const F_IntersectPoint * & ip = _grid->_gridIntP[ hex->_origNodeInd + _nodeShift[iN] ];
3608           if ( !ip )
3609           {
3610             ip = _grid->_extIntPool.getNew();
3611             ip->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
3612             //ip->_transition = Trans_INTERNAL;
3613           }
3614           else if ( ip->_faceIDs.back() != _grid->PseudoIntExtFaceID() )
3615           {
3616             ip->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
3617           }
3618           hex->_nbFaceIntNodes++;
3619           //isCut = true;
3620         }
3621         else
3622         {
3623           const gp_Pnt&      p = intersection.Point( 1 );
3624           F_IntersectPoint* ip = _grid->_extIntPool.getNew();
3625           ip->_node = meshDS->AddNode( p.X(), p.Y(), p.Z() );
3626           ip->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
3627           ip->_transition = Trans_INTERNAL;
3628           meshDS->SetNodeInVolume( ip->_node, solidID );
3629
3630           CellsAroundLink fourCells( _grid, iDir );
3631           fourCells.Init( hex->_i, hex->_j, hex->_k, iLink );
3632           int i,j,k, cellIndex;
3633           for ( int iC = 0; iC < 4; ++iC ) // loop on 4 cells sharing the link
3634           {
3635             if ( !fourCells.GetCell( iC, i,j,k, cellIndex ))
3636               continue;
3637             Hexahedron * h = hexes[ cellIndex ];
3638             if ( !h )
3639               h = hexes[ cellIndex ] = new Hexahedron( *this, i, j, k, cellIndex );
3640             const int iL = iC + iDir * 4;
3641             h->_hexLinks[iL]._fIntPoints.push_back( ip );
3642             h->_nbFaceIntNodes++;
3643             //isCut = true;
3644           }
3645         }
3646       }
3647
3648       // if ( isCut )
3649       //   for ( size_t i = 0; i < hex->_eIntPoints.size(); ++i )
3650       //   {
3651       //     if ( _grid->IsInternal( hex->_eIntPoints[i]->_shapeID ) &&
3652       //          ! hex->_eIntPoints[i]->IsOnFace( _grid->PseudoIntExtFaceID() ))
3653       //       hex->_eIntPoints[i]->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
3654       //   }
3655       continue;
3656
3657     } // loop on all hexes
3658     return;
3659   }
3660
3661   //================================================================================
3662   /*!
3663    * \brief Return intersection point on INTERNAL FACE most distant from given ones
3664    */
3665   gp_Pnt Hexahedron::mostDistantInternalPnt( int hexIndex, const gp_Pnt& p1, const gp_Pnt& p2 )
3666   {
3667     gp_Pnt resultPnt = p1;
3668
3669     double maxDist2 = 0;
3670     for ( int iLink = 0; iLink < 12; ++iLink ) // check links
3671     {
3672       _Link& link = _hexLinks[ iLink ];
3673       for ( size_t i = 0; i < link._fIntPoints.size(); ++i )
3674         if ( _grid->PseudoIntExtFaceID() != link._fIntPoints[i]->_faceIDs[0] &&
3675              _grid->IsInternal( link._fIntPoints[i]->_faceIDs[0] ) &&
3676              link._fIntPoints[i]->_node )
3677         {
3678           gp_Pnt p = SMESH_NodeXYZ( link._fIntPoints[i]->_node );
3679           double d = p1.SquareDistance( p );
3680           if ( d > maxDist2 )
3681           {
3682             resultPnt = p;
3683             maxDist2  = d;
3684           }
3685           else
3686           {
3687             d = p2.SquareDistance( p );
3688             if ( d > maxDist2 )
3689             {
3690               resultPnt = p;
3691               maxDist2  = d;
3692             }
3693           }
3694         }
3695     }
3696     setIJK( hexIndex );
3697     _origNodeInd = _grid->NodeIndex( _i,_j,_k );
3698
3699     for ( size_t iN = 0; iN < 8; ++iN ) // check corners
3700     {
3701       _hexNodes[iN]._node     = _grid->_nodes   [ _origNodeInd + _nodeShift[iN] ];
3702       _hexNodes[iN]._intPoint = _grid->_gridIntP[ _origNodeInd + _nodeShift[iN] ];
3703       if ( _hexNodes[iN]._intPoint )
3704         for ( size_t iF = 0; iF < _hexNodes[iN]._intPoint->_faceIDs.size(); ++iF )
3705         {
3706           if ( _grid->IsInternal( _hexNodes[iN]._intPoint->_faceIDs[iF]))
3707           {
3708             gp_Pnt p = SMESH_NodeXYZ( _hexNodes[iN]._node );
3709             double d = p1.SquareDistance( p );
3710             if ( d > maxDist2 )
3711             {
3712               resultPnt = p;
3713               maxDist2  = d;
3714             }
3715             else
3716             {
3717               d = p2.SquareDistance( p );
3718               if ( d > maxDist2 )
3719               {
3720                 resultPnt = p;
3721                 maxDist2  = d;
3722               }
3723             }
3724           }
3725         }
3726     }
3727     if ( maxDist2 < _grid->_tol * _grid->_tol )
3728       return p1;
3729
3730     return resultPnt;
3731   }
3732
3733   //================================================================================
3734   /*!
3735    * \brief Finds intersection of a curve with a plane
3736    *  \param [in] u1 - parameter of one curve point
3737    *  \param [in] proj1 - projection of the curve point to the plane normal
3738    *  \param [in] u2 - parameter of another curve point
3739    *  \param [in] proj2 - projection of the other curve point to the plane normal
3740    *  \param [in] proj - projection of a point where the curve intersects the plane
3741    *  \param [in] curve - the curve
3742    *  \param [in] axis - the plane normal
3743    *  \param [in] origin - the plane origin
3744    *  \return gp_Pnt - the found intersection point
3745    */
3746   gp_Pnt Hexahedron::findIntPoint( double u1, double proj1,
3747                                    double u2, double proj2,
3748                                    double proj,
3749                                    BRepAdaptor_Curve& curve,
3750                                    const gp_XYZ& axis,
3751                                    const gp_XYZ& origin)
3752   {
3753     double r = (( proj - proj1 ) / ( proj2 - proj1 ));
3754     double u = u1 * ( 1 - r ) + u2 * r;
3755     gp_Pnt p = curve.Value( u );
3756     double newProj =  axis * ( p.XYZ() - origin );
3757     if ( Abs( proj - newProj ) > _grid->_tol / 10. )
3758     {
3759       if ( r > 0.5 )
3760         return findIntPoint( u2, proj2, u, newProj, proj, curve, axis, origin );
3761       else
3762         return findIntPoint( u1, proj2, u, newProj, proj, curve, axis, origin );
3763     }
3764     return p;
3765   }
3766
3767   //================================================================================
3768   /*!
3769    * \brief Returns indices of a hexahedron sub-entities holding a point
3770    *  \param [in] ip - intersection point
3771    *  \param [out] facets - 0-3 facets holding a point
3772    *  \param [out] sub - index of a vertex or an edge holding a point
3773    *  \return int - number of facets holding a point
3774    */
3775   int Hexahedron::getEntity( const E_IntersectPoint* ip, int* facets, int& sub )
3776   {
3777     enum { X = 1, Y = 2, Z = 4 }; // == 001, 010, 100
3778     int nbFacets = 0;
3779     int vertex = 0, edgeMask = 0;
3780
3781     if ( Abs( _grid->_coords[0][ _i   ] - ip->_uvw[0] ) < _grid->_tol ) {
3782       facets[ nbFacets++ ] = SMESH_Block::ID_F0yz;
3783       edgeMask |= X;
3784     }
3785     else if ( Abs( _grid->_coords[0][ _i+1 ] - ip->_uvw[0] ) < _grid->_tol ) {
3786       facets[ nbFacets++ ] = SMESH_Block::ID_F1yz;
3787       vertex   |= X;
3788       edgeMask |= X;
3789     }
3790     if ( Abs( _grid->_coords[1][ _j   ] - ip->_uvw[1] ) < _grid->_tol ) {
3791       facets[ nbFacets++ ] = SMESH_Block::ID_Fx0z;
3792       edgeMask |= Y;
3793     }
3794     else if ( Abs( _grid->_coords[1][ _j+1 ] - ip->_uvw[1] ) < _grid->_tol ) {
3795       facets[ nbFacets++ ] = SMESH_Block::ID_Fx1z;
3796       vertex   |= Y;
3797       edgeMask |= Y;
3798     }
3799     if ( Abs( _grid->_coords[2][ _k   ] - ip->_uvw[2] ) < _grid->_tol ) {
3800       facets[ nbFacets++ ] = SMESH_Block::ID_Fxy0;
3801       edgeMask |= Z;
3802     }
3803     else if ( Abs( _grid->_coords[2][ _k+1 ] - ip->_uvw[2] ) < _grid->_tol ) {
3804       facets[ nbFacets++ ] = SMESH_Block::ID_Fxy1;
3805       vertex   |= Z;
3806       edgeMask |= Z;
3807     }
3808
3809     switch ( nbFacets )
3810     {
3811     case 0: sub = 0;         break;
3812     case 1: sub = facets[0]; break;
3813     case 2: {
3814       const int edge [3][8] = {
3815         { SMESH_Block::ID_E00z, SMESH_Block::ID_E10z,
3816           SMESH_Block::ID_E01z, SMESH_Block::ID_E11z },
3817         { SMESH_Block::ID_E0y0, SMESH_Block::ID_E1y0, 0, 0,
3818           SMESH_Block::ID_E0y1, SMESH_Block::ID_E1y1 },
3819         { SMESH_Block::ID_Ex00, 0, SMESH_Block::ID_Ex10, 0,
3820           SMESH_Block::ID_Ex01, 0, SMESH_Block::ID_Ex11 }
3821       };
3822       switch ( edgeMask ) {
3823       case X | Y: sub = edge[ 0 ][ vertex ]; break;
3824       case X | Z: sub = edge[ 1 ][ vertex ]; break;
3825       default:    sub = edge[ 2 ][ vertex ];
3826       }
3827       break;
3828     }
3829     //case 3:
3830     default:
3831       sub = vertex + SMESH_Block::ID_FirstV;
3832     }
3833
3834     return nbFacets;
3835   }
3836   //================================================================================
3837   /*!
3838    * \brief Adds intersection with an EDGE
3839    */
3840   bool Hexahedron::addIntersection( const E_IntersectPoint* ip,
3841                                     vector< Hexahedron* >&  hexes,
3842                                     int ijk[], int dIJK[] )
3843   {
3844     bool added = false;
3845
3846     size_t hexIndex[4] = {
3847       _grid->CellIndex( ijk[0], ijk[1], ijk[2] ),
3848       dIJK[0] ? _grid->CellIndex( ijk[0]+dIJK[0], ijk[1], ijk[2] ) : -1,
3849       dIJK[1] ? _grid->CellIndex( ijk[0], ijk[1]+dIJK[1], ijk[2] ) : -1,
3850       dIJK[2] ? _grid->CellIndex( ijk[0], ijk[1], ijk[2]+dIJK[2] ) : -1
3851     };
3852     for ( int i = 0; i < 4; ++i )
3853     {
3854       if ( hexIndex[i] < hexes.size() && hexes[ hexIndex[i] ] )
3855       {
3856         Hexahedron* h = hexes[ hexIndex[i] ];
3857         h->_eIntPoints.reserve(2);
3858         h->_eIntPoints.push_back( ip );
3859         added = true;
3860 #ifdef _DEBUG_
3861         // check if ip is really inside the hex
3862         if ( h->isOutParam( ip->_uvw ))
3863           throw SALOME_Exception("ip outside a hex");
3864 #endif
3865       }
3866     }
3867     return added;
3868   }
3869   //================================================================================
3870   /*!
3871    * \brief Finds nodes at a path from one node to another via intersections with EDGEs
3872    */
3873   bool Hexahedron::findChain( _Node*          n1,
3874                               _Node*          n2,
3875                               _Face&          quad,
3876                               vector<_Node*>& chn )
3877   {
3878     chn.clear();
3879     chn.push_back( n1 );
3880     for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
3881       if ( !quad._eIntNodes[ iP ]->IsUsedInFace( &quad ) &&
3882            n1->IsLinked( quad._eIntNodes[ iP ]->_intPoint ) &&
3883            n2->IsLinked( quad._eIntNodes[ iP ]->_intPoint ))
3884       {
3885         chn.push_back( quad._eIntNodes[ iP ]);
3886         chn.push_back( n2 );
3887         quad._eIntNodes[ iP ]->_usedInFace = &quad;
3888         return true;
3889       }
3890     bool found;
3891     do
3892     {
3893       found = false;
3894       for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
3895         if ( !quad._eIntNodes[ iP ]->IsUsedInFace( &quad ) &&
3896              chn.back()->IsLinked( quad._eIntNodes[ iP ]->_intPoint ))
3897         {
3898           chn.push_back( quad._eIntNodes[ iP ]);
3899           found = ( quad._eIntNodes[ iP ]->_usedInFace = &quad );
3900           break;
3901         }
3902     } while ( found && ! chn.back()->IsLinked( n2->_intPoint ) );
3903
3904     if ( chn.back() != n2 && chn.back()->IsLinked( n2->_intPoint ))
3905       chn.push_back( n2 );
3906
3907     return chn.size() > 1;
3908   }
3909   //================================================================================
3910   /*!
3911    * \brief Try to heal a polygon whose ends are not connected
3912    */
3913   bool Hexahedron::closePolygon( _Face* polygon, vector<_Node*>& chainNodes ) const
3914   {
3915     int i = -1, nbLinks = polygon->_links.size();
3916     if ( nbLinks < 3 )
3917       return false;
3918     vector< _OrientedLink > newLinks;
3919     // find a node lying on the same FACE as the last one
3920     _Node*   node = polygon->_links.back().LastNode();
3921     int avoidFace = node->IsLinked( polygon->_links.back().FirstNode()->_intPoint );
3922     for ( i = nbLinks - 2; i >= 0; --i )
3923       if ( node->IsLinked( polygon->_links[i].FirstNode()->_intPoint, avoidFace ))
3924         break;
3925     if ( i >= 0 )
3926     {
3927       for ( ; i < nbLinks; ++i )
3928         newLinks.push_back( polygon->_links[i] );
3929     }
3930     else
3931     {
3932       // find a node lying on the same FACE as the first one
3933       node      = polygon->_links[0].FirstNode();
3934       avoidFace = node->IsLinked( polygon->_links[0].LastNode()->_intPoint );
3935       for ( i = 1; i < nbLinks; ++i )
3936         if ( node->IsLinked( polygon->_links[i].LastNode()->_intPoint, avoidFace ))
3937           break;
3938       if ( i < nbLinks )
3939         for ( nbLinks = i + 1, i = 0; i < nbLinks; ++i )
3940           newLinks.push_back( polygon->_links[i] );
3941     }
3942     if ( newLinks.size() > 1 )
3943     {
3944       polygon->_links.swap( newLinks );
3945       chainNodes.clear();
3946       chainNodes.push_back( polygon->_links.back().LastNode() );
3947       chainNodes.push_back( polygon->_links[0].FirstNode() );
3948       return true;
3949     }
3950     return false;
3951   }
3952   //================================================================================
3953   /*!
3954    * \brief Finds nodes on the same EDGE as the first node of avoidSplit.
3955    *
3956    * This function is for
3957    * 1) a case where an EDGE lies on a quad which lies on a FACE
3958    *    so that a part of quad in ON and another part is IN
3959    * 2) INTERNAL FACE passes through the 1st node of avoidSplit
3960    */
3961   bool Hexahedron::findChainOnEdge( const vector< _OrientedLink >& splits,
3962                                     const _OrientedLink&           prevSplit,
3963                                     const _OrientedLink&           avoidSplit,
3964                                     size_t &                       iS,
3965                                     _Face&                         quad,
3966                                     vector<_Node*>&                chn )
3967   {
3968     _Node* pn1 = prevSplit.FirstNode();
3969     _Node* pn2 = prevSplit.LastNode();
3970     int avoidFace = pn1->IsLinked( pn2->_intPoint ); // FACE under the quad
3971     if ( avoidFace < 1 && pn1->_intPoint )
3972       return false;
3973
3974     _Node* n = 0, *stopNode = avoidSplit.LastNode();
3975
3976     chn.clear();
3977     if ( !quad._eIntNodes.empty() ) // connect pn2 with EDGE intersections
3978     {
3979       chn.push_back( pn2 );
3980       bool found;
3981       do
3982       {
3983         found = false;
3984         for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
3985           if (( !quad._eIntNodes[ iP ]->IsUsedInFace( &quad )) &&
3986               ( chn.back()->IsLinked( quad._eIntNodes[ iP ]->_intPoint, avoidFace )) &&
3987               ( !avoidFace || quad._eIntNodes[ iP ]->IsOnFace( avoidFace )))
3988           {
3989             chn.push_back( quad._eIntNodes[ iP ]);
3990             found = ( quad._eIntNodes[ iP ]->_usedInFace = &quad );
3991             break;
3992           }
3993       } while ( found );
3994       pn2 = chn.back();
3995     }
3996
3997     int i;
3998     for ( i = splits.size()-1; i >= 0; --i ) // connect new pn2 (at _eIntNodes) with a split
3999     {
4000       if ( !splits[i] )
4001         continue;
4002
4003       n = splits[i].LastNode();
4004       if ( n == stopNode )
4005         break;
4006       if (( n != pn1 ) &&
4007           ( n->IsLinked( pn2->_intPoint, avoidFace )) &&
4008           ( !avoidFace || n->IsOnFace( avoidFace )))
4009         break;
4010
4011       n = splits[i].FirstNode();
4012       if ( n == stopNode )
4013         break;
4014       if (( n->IsLinked( pn2->_intPoint, avoidFace )) &&
4015           ( !avoidFace || n->IsOnFace( avoidFace )))
4016         break;
4017       n = 0;
4018     }
4019     if ( n && n != stopNode )
4020     {
4021       if ( chn.empty() )
4022         chn.push_back( pn2 );
4023       chn.push_back( n );
4024       iS = i-1;
4025       return true;
4026     }
4027     else if ( !chn.empty() && chn.back()->_isInternalFlags )
4028     {
4029       // INTERNAL FACE partially cuts the quad
4030       for ( int i = chn.size() - 2; i >= 0; --i )
4031         chn.push_back( chn[ i ]);
4032       return true;
4033     }
4034     return false;
4035   }
4036   //================================================================================
4037   /*!
4038    * \brief Checks transition at the ginen intersection node of a link
4039    */
4040   bool Hexahedron::isOutPoint( _Link& link, int iP,
4041                                SMESH_MesherHelper& helper, const Solid* solid ) const
4042   {
4043     bool isOut = false;
4044
4045     if ( link._fIntNodes[iP]->faces().size() == 1 &&
4046          _grid->IsInternal( link._fIntNodes[iP]->face(0) ))
4047       return false;
4048
4049     const bool moreIntPoints = ( iP+1 < (int) link._fIntNodes.size() );
4050
4051     // get 2 _Node's
4052     _Node* n1 = link._fIntNodes[ iP ];
4053     if ( !n1->Node() )
4054       n1 = link._nodes[0];
4055     _Node* n2 = moreIntPoints ? link._fIntNodes[ iP+1 ] : 0;
4056     if ( !n2 || !n2->Node() )
4057       n2 = link._nodes[1];
4058     if ( !n2->Node() )
4059       return true;
4060
4061     // get all FACEs under n1 and n2
4062     set< TGeomID > faceIDs;
4063     if ( moreIntPoints ) faceIDs.insert( link._fIntNodes[iP+1]->faces().begin(),
4064                                          link._fIntNodes[iP+1]->faces().end() );
4065     if ( n2->_intPoint ) faceIDs.insert( n2->_intPoint->_faceIDs.begin(),
4066                                          n2->_intPoint->_faceIDs.end() );
4067     if ( faceIDs.empty() )
4068       return false; // n2 is inside
4069     if ( n1->_intPoint ) faceIDs.insert( n1->_intPoint->_faceIDs.begin(),
4070                                          n1->_intPoint->_faceIDs.end() );
4071     faceIDs.insert( link._fIntNodes[iP]->faces().begin(),
4072                     link._fIntNodes[iP]->faces().end() );
4073
4074     // get a point between 2 nodes
4075     gp_Pnt p1      = n1->Point();
4076     gp_Pnt p2      = n2->Point();
4077     gp_Pnt pOnLink = 0.8 * p1.XYZ() + 0.2 * p2.XYZ();
4078
4079     TopLoc_Location loc;
4080
4081     set< TGeomID >::iterator faceID = faceIDs.begin();
4082     for ( ; faceID != faceIDs.end(); ++faceID )
4083     {
4084       // project pOnLink on a FACE
4085       if ( *faceID < 1 || !solid->Contains( *faceID )) continue;
4086       const TopoDS_Face& face = TopoDS::Face( _grid->Shape( *faceID ));
4087       GeomAPI_ProjectPointOnSurf& proj = helper.GetProjector( face, loc, 0.1*_grid->_tol );
4088       gp_Pnt testPnt = pOnLink.Transformed( loc.Transformation().Inverted() );
4089       proj.Perform( testPnt );
4090       if ( proj.IsDone() && proj.NbPoints() > 0 )       
4091       {
4092         Standard_Real u,v;
4093         proj.LowerDistanceParameters( u,v );
4094
4095         if ( proj.LowerDistance() <= 0.1 * _grid->_tol )
4096         {
4097           isOut = false;
4098         }
4099         else
4100         {
4101           // find isOut by normals
4102           gp_Dir normal;
4103           if ( GeomLib::NormEstim( BRep_Tool::Surface( face, loc ),
4104                                    gp_Pnt2d( u,v ),
4105                                    0.1*_grid->_tol,
4106                                    normal ) < 3 )
4107           {
4108             if ( solid->Orientation( face ) == TopAbs_REVERSED )
4109               normal.Reverse();
4110             gp_Vec v( proj.NearestPoint(), testPnt );
4111             isOut = ( v * normal > 0 );
4112           }
4113         }
4114         if ( !isOut )
4115         {
4116           // classify a projection
4117           if ( !n1->IsOnFace( *faceID ) || !n2->IsOnFace( *faceID ))
4118           {
4119             BRepTopAdaptor_FClass2d cls( face, Precision::Confusion() );
4120             TopAbs_State state = cls.Perform( gp_Pnt2d( u,v ));
4121             if ( state == TopAbs_OUT )
4122             {
4123               isOut = true;
4124               continue;
4125             }
4126           }
4127           return false;
4128         }
4129       }
4130     }
4131     return isOut;
4132   }
4133   //================================================================================
4134   /*!
4135    * \brief Sort nodes on a FACE
4136    */
4137   void Hexahedron::sortVertexNodes(vector<_Node*>& nodes, _Node* curNode, TGeomID faceID)
4138   {
4139     if ( nodes.size() > 20 ) return;
4140
4141     // get shapes under nodes
4142     TGeomID nShapeIds[20], *nShapeIdsEnd = &nShapeIds[0] + nodes.size();
4143     for ( size_t i = 0; i < nodes.size(); ++i )
4144       if ( !( nShapeIds[i] = nodes[i]->ShapeID() ))
4145         return;
4146
4147     // get shapes of the FACE
4148     const TopoDS_Face&  face = TopoDS::Face( _grid->Shape( faceID ));
4149     list< TopoDS_Edge > edges;
4150     list< int >         nbEdges;
4151     int nbW = SMESH_Block::GetOrderedEdges (face, edges, nbEdges);
4152     if ( nbW > 1 ) {
4153       // select a WIRE - remove EDGEs of irrelevant WIREs from edges
4154       list< TopoDS_Edge >::iterator e = edges.begin(), eEnd = e;
4155       list< int >::iterator nE = nbEdges.begin();
4156       for ( ; nbW > 0; ++nE, --nbW )
4157       {
4158         std::advance( eEnd, *nE );
4159         for ( ; e != eEnd; ++e )
4160           for ( int i = 0; i < 2; ++i )
4161           {
4162             TGeomID id = i==0 ?
4163               _grid->ShapeID( *e ) :
4164               _grid->ShapeID( SMESH_MesherHelper::IthVertex( 0, *e ));
4165             if (( id > 0 ) &&
4166                 ( std::find( &nShapeIds[0], nShapeIdsEnd, id ) != nShapeIdsEnd ))
4167             {
4168               edges.erase( eEnd, edges.end() ); // remove rest wires
4169               e = eEnd = edges.end();
4170               --e;
4171               nbW = 0;
4172               break;
4173             }
4174           }
4175         if ( nbW > 0 )
4176           edges.erase( edges.begin(), eEnd ); // remove a current irrelevant wire
4177       }
4178     }
4179     // rotate edges to have the first one at least partially out of the hexa
4180     list< TopoDS_Edge >::iterator e = edges.begin(), eMidOut = edges.end();
4181     for ( ; e != edges.end(); ++e )
4182     {
4183       if ( !_grid->ShapeID( *e ))
4184         continue;
4185       bool isOut = false;
4186       gp_Pnt p;
4187       double uvw[3], f,l;
4188       for ( int i = 0; i < 2 && !isOut; ++i )
4189       {
4190         if ( i == 0 )
4191         {
4192           TopoDS_Vertex v = SMESH_MesherHelper::IthVertex( 0, *e );
4193           p = BRep_Tool::Pnt( v );
4194         }
4195         else if ( eMidOut == edges.end() )
4196         {
4197           TopLoc_Location loc;
4198           Handle(Geom_Curve) c = BRep_Tool::Curve( *e, loc, f, l);
4199           if ( c.IsNull() ) break;
4200           p = c->Value( 0.5 * ( f + l )).Transformed( loc );
4201         }
4202         else
4203         {
4204           continue;
4205         }
4206
4207         _grid->ComputeUVW( p.XYZ(), uvw );
4208         if ( isOutParam( uvw ))
4209         {
4210           if ( i == 0 )
4211             isOut = true;
4212           else
4213             eMidOut = e;
4214         }
4215       }
4216       if ( isOut )
4217         break;
4218     }
4219     if ( e != edges.end() )
4220       edges.splice( edges.end(), edges, edges.begin(), e );
4221     else if ( eMidOut != edges.end() )
4222       edges.splice( edges.end(), edges, edges.begin(), eMidOut );
4223
4224     // sort nodes according to the order of edges
4225     _Node*  orderNodes   [20];
4226     //TGeomID orderShapeIDs[20];
4227     size_t nbN = 0;
4228     TGeomID id, *pID = 0;
4229     for ( e = edges.begin(); e != edges.end(); ++e )
4230     {
4231       if (( id = _grid->ShapeID( SMESH_MesherHelper::IthVertex( 0, *e ))) &&
4232           (( pID = std::find( &nShapeIds[0], nShapeIdsEnd, id )) != nShapeIdsEnd ))
4233       {
4234         //orderShapeIDs[ nbN ] = id;
4235         orderNodes   [ nbN++ ] = nodes[ pID - &nShapeIds[0] ];
4236         *pID = -1;
4237       }
4238       if (( id = _grid->ShapeID( *e )) &&
4239           (( pID = std::find( &nShapeIds[0], nShapeIdsEnd, id )) != nShapeIdsEnd ))
4240       {
4241         //orderShapeIDs[ nbN ] = id;
4242         orderNodes   [ nbN++ ] = nodes[ pID - &nShapeIds[0] ];
4243         *pID = -1;
4244       }
4245     }
4246     if ( nbN != nodes.size() )
4247       return;
4248
4249     bool reverse = ( orderNodes[0    ]->Point().SquareDistance( curNode->Point() ) >
4250                      orderNodes[nbN-1]->Point().SquareDistance( curNode->Point() ));
4251
4252     for ( size_t i = 0; i < nodes.size(); ++i )
4253       nodes[ i ] = orderNodes[ reverse ? nbN-1-i : i ];
4254   }
4255
4256   //================================================================================
4257   /*!
4258    * \brief Adds computed elements to the mesh
4259    */
4260   int Hexahedron::addVolumes( SMESH_MesherHelper& helper )
4261   {
4262     F_IntersectPoint noIntPnt;
4263     const bool toCheckNodePos = _grid->IsToCheckNodePos();
4264
4265     int nbAdded = 0;
4266     // add elements resulted from hexahedron intersection
4267     for ( _volumeDef* volDef = &_volumeDefs; volDef; volDef = volDef->_next )
4268     {
4269       vector< const SMDS_MeshNode* > nodes( volDef->_nodes.size() );
4270       for ( size_t iN = 0; iN < nodes.size(); ++iN )
4271       {
4272         if ( !( nodes[iN] = volDef->_nodes[iN].Node() ))
4273         {
4274           if ( const E_IntersectPoint* eip = volDef->_nodes[iN].EdgeIntPnt() )
4275           {
4276             nodes[iN] = volDef->_nodes[iN]._intPoint->_node =
4277               helper.AddNode( eip->_point.X(),
4278                               eip->_point.Y(),
4279                               eip->_point.Z() );
4280             if ( _grid->ShapeType( eip->_shapeID ) == TopAbs_VERTEX )
4281               helper.GetMeshDS()->SetNodeOnVertex( nodes[iN], eip->_shapeID );
4282             else
4283               helper.GetMeshDS()->SetNodeOnEdge( nodes[iN], eip->_shapeID );
4284           }
4285           else
4286             throw SALOME_Exception("Bug: no node at intersection point");
4287         }
4288         else if ( volDef->_nodes[iN]._intPoint &&
4289                   volDef->_nodes[iN]._intPoint->_node == volDef->_nodes[iN]._node )
4290         {
4291           // Update position of node at EDGE intersection;
4292           // see comment to _Node::Add( E_IntersectPoint )
4293           SMESHDS_Mesh* mesh = helper.GetMeshDS();
4294           TGeomID    shapeID = volDef->_nodes[iN].EdgeIntPnt()->_shapeID;
4295           mesh->UnSetNodeOnShape( nodes[iN] );
4296           if ( _grid->ShapeType( shapeID ) == TopAbs_VERTEX )
4297             mesh->SetNodeOnVertex( nodes[iN], shapeID );
4298           else
4299             mesh->SetNodeOnEdge( nodes[iN], shapeID );
4300         }
4301         else if ( toCheckNodePos &&
4302                   !nodes[iN]->isMarked() && 
4303                   _grid->ShapeType( nodes[iN]->GetShapeID() ) == TopAbs_FACE )
4304         {
4305           _grid->SetOnShape( nodes[iN], noIntPnt, /*unset=*/true );
4306           nodes[iN]->setIsMarked( true );
4307         }
4308       }
4309
4310       const SMDS_MeshElement* v = 0;
4311       if ( !volDef->_quantities.empty() )
4312       {
4313         v = helper.AddPolyhedralVolume( nodes, volDef->_quantities );
4314       }
4315       else
4316       {
4317         switch ( nodes.size() )
4318         {
4319         case 8: v = helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],
4320                                       nodes[4],nodes[5],nodes[6],nodes[7] );
4321           break;
4322         case 4: v = helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3] );
4323           break;
4324         case 6: v = helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],nodes[4],nodes[5] );
4325           break;
4326         case 5: v = helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],nodes[4] );
4327           break;
4328         }
4329       }
4330       if (( volDef->_volume = v ))
4331       {
4332         helper.GetMeshDS()->SetMeshElementOnShape( v, volDef->_solidID );
4333         ++nbAdded;
4334       }
4335     }
4336
4337     return nbAdded;
4338   }
4339   //================================================================================
4340   /*!
4341    * \brief Return true if the element is in a hole
4342    */
4343   bool Hexahedron::isInHole() const
4344   {
4345     if ( !_vIntNodes.empty() )
4346       return false;
4347
4348     const size_t ijk[3] = { _i, _j, _k };
4349     F_IntersectPoint curIntPnt;
4350
4351     // consider a cell to be in a hole if all links in any direction
4352     // comes OUT of geometry
4353     for ( int iDir = 0; iDir < 3; ++iDir )
4354     {
4355       const vector<double>& coords = _grid->_coords[ iDir ];
4356       LineIndexer               li = _grid->GetLineIndexer( iDir );
4357       li.SetIJK( _i,_j,_k );
4358       size_t lineIndex[4] = { li.LineIndex  (),
4359                               li.LineIndex10(),
4360                               li.LineIndex01(),
4361                               li.LineIndex11() };
4362       bool allLinksOut = true, hasLinks = false;
4363       for ( int iL = 0; iL < 4 && allLinksOut; ++iL ) // loop on 4 links parallel to iDir
4364       {
4365         const _Link& link = _hexLinks[ iL + 4*iDir ];
4366         // check transition of the first node of a link
4367         const F_IntersectPoint* firstIntPnt = 0;
4368         if ( link._nodes[0]->Node() ) // 1st node is a hexa corner
4369         {
4370           curIntPnt._paramOnLine = coords[ ijk[ iDir ]] - coords[0] + _grid->_tol;
4371           const GridLine& line = _grid->_lines[ iDir ][ lineIndex[ iL ]];
4372           multiset< F_IntersectPoint >::const_iterator ip =
4373             line._intPoints.upper_bound( curIntPnt );
4374           --ip;
4375           firstIntPnt = &(*ip);
4376         }
4377         else if ( !link._fIntPoints.empty() )
4378         {
4379           firstIntPnt = link._fIntPoints[0];
4380         }
4381
4382         if ( firstIntPnt )
4383         {
4384           hasLinks = true;
4385           allLinksOut = ( firstIntPnt->_transition == Trans_OUT &&
4386                           !_grid->IsShared( firstIntPnt->_faceIDs[0] ));
4387         }
4388       }
4389       if ( hasLinks && allLinksOut )
4390         return true;
4391     }
4392     return false;
4393   }
4394
4395   //================================================================================
4396   /*!
4397    * \brief Check if a polyherdon has an edge lying on EDGE shared by strange FACE
4398    *        that will be meshed by other algo
4399    */
4400   bool Hexahedron::hasStrangeEdge() const
4401   {
4402     if ( _eIntPoints.size() < 2 )
4403       return false;
4404
4405     TopTools_MapOfShape edges;
4406     for ( size_t i = 0; i < _eIntPoints.size(); ++i )
4407     {
4408       if ( !_grid->IsStrangeEdge( _eIntPoints[i]->_shapeID ))
4409         continue;
4410       const TopoDS_Shape& s = _grid->Shape( _eIntPoints[i]->_shapeID );
4411       if ( s.ShapeType() == TopAbs_EDGE )
4412       {
4413         if ( ! edges.Add( s ))
4414           return true; // an EDGE encounters twice
4415       }
4416       else
4417       {
4418         PShapeIteratorPtr edgeIt = _grid->_helper->GetAncestors( s,
4419                                                                  *_grid->_helper->GetMesh(),
4420                                                                  TopAbs_EDGE );
4421         while ( const TopoDS_Shape* edge = edgeIt->next() )
4422           if ( ! edges.Add( *edge ))
4423             return true; // an EDGE encounters twice
4424       }
4425     }
4426     return false;
4427   }
4428
4429   //================================================================================
4430   /*!
4431    * \brief Return true if a polyhedron passes _sizeThreshold criterion
4432    */
4433   bool Hexahedron::checkPolyhedronSize( bool cutByInternalFace ) const
4434   {
4435     if ( cutByInternalFace && !_grid->_toUseThresholdForInternalFaces )
4436     {
4437       // check if any polygon fully lies on shared/internal FACEs
4438       for ( size_t iP = 0; iP < _polygons.size(); ++iP )
4439       {
4440         const _Face& polygon = _polygons[iP];
4441         if ( polygon._links.empty() )
4442           continue;
4443         bool allNodesInternal = true;
4444         for ( size_t iL = 0; iL < polygon._links.size() &&  allNodesInternal; ++iL )
4445         {
4446           _Node* n = polygon._links[ iL ].FirstNode();
4447           allNodesInternal = (( n->IsCutByInternal() ) ||
4448                               ( n->_intPoint && _grid->IsAnyShared( n->_intPoint->_faceIDs )));
4449         }
4450         if ( allNodesInternal )
4451           return true;
4452       }
4453     }
4454     if ( this->hasStrangeEdge() )
4455       return true;
4456
4457     double volume = 0;
4458     for ( size_t iP = 0; iP < _polygons.size(); ++iP )
4459     {
4460       const _Face& polygon = _polygons[iP];
4461       if ( polygon._links.empty() )
4462         continue;
4463       gp_XYZ area (0,0,0);
4464       gp_XYZ p1 = polygon._links[ 0 ].FirstNode()->Point().XYZ();
4465       for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
4466       {
4467         gp_XYZ p2 = polygon._links[ iL ].LastNode()->Point().XYZ();
4468         area += p1 ^ p2;
4469         p1 = p2;
4470       }
4471       volume += p1 * area;
4472     }
4473     volume /= 6;
4474
4475     double initVolume = _sideLength[0] * _sideLength[1] * _sideLength[2];
4476
4477     return volume > initVolume / _grid->_sizeThreshold;
4478   }
4479   //================================================================================
4480   /*!
4481    * \brief Tries to create a hexahedron
4482    */
4483   bool Hexahedron::addHexa()
4484   {
4485     int nbQuad = 0, iQuad = -1;
4486     for ( size_t i = 0; i < _polygons.size(); ++i )
4487     {
4488       if ( _polygons[i]._links.empty() )
4489         continue;
4490       if ( _polygons[i]._links.size() != 4 )
4491         return false;
4492       ++nbQuad;
4493       if ( iQuad < 0 )
4494         iQuad = i;
4495     }
4496     if ( nbQuad != 6 )
4497       return false;
4498
4499     _Node* nodes[8];
4500     int nbN = 0;
4501     for ( int iL = 0; iL < 4; ++iL )
4502     {
4503       // a base node
4504       nodes[iL] = _polygons[iQuad]._links[iL].FirstNode();
4505       ++nbN;
4506
4507       // find a top node above the base node
4508       _Link* link = _polygons[iQuad]._links[iL]._link;
4509       if ( !link->_faces[0] || !link->_faces[1] )
4510         return debugDumpLink( link );
4511       // a quadrangle sharing <link> with _polygons[iQuad]
4512       _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[iQuad] )];
4513       for ( int i = 0; i < 4; ++i )
4514         if ( quad->_links[i]._link == link )
4515         {
4516           // 1st node of a link opposite to <link> in <quad>
4517           nodes[iL+4] = quad->_links[(i+2)%4].FirstNode();
4518           ++nbN;
4519           break;
4520         }
4521     }
4522     if ( nbN == 8 )
4523       _volumeDefs.Set( &nodes[0], 8 );
4524
4525     return nbN == 8;
4526   }
4527   //================================================================================
4528   /*!
4529    * \brief Tries to create a tetrahedron
4530    */
4531   bool Hexahedron::addTetra()
4532   {
4533     int iTria = -1;
4534     for ( size_t i = 0; i < _polygons.size() && iTria < 0; ++i )
4535       if ( _polygons[i]._links.size() == 3 )
4536         iTria = i;
4537     if ( iTria < 0 )
4538       return false;
4539
4540     _Node* nodes[4];
4541     nodes[0] = _polygons[iTria]._links[0].FirstNode();
4542     nodes[1] = _polygons[iTria]._links[1].FirstNode();
4543     nodes[2] = _polygons[iTria]._links[2].FirstNode();
4544
4545     _Link* link = _polygons[iTria]._links[0]._link;
4546     if ( !link->_faces[0] || !link->_faces[1] )
4547       return debugDumpLink( link );
4548
4549     // a triangle sharing <link> with _polygons[0]
4550     _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[iTria] )];
4551     for ( int i = 0; i < 3; ++i )
4552       if ( tria->_links[i]._link == link )
4553       {
4554         nodes[3] = tria->_links[(i+1)%3].LastNode();
4555         _volumeDefs.Set( &nodes[0], 4 );
4556         return true;
4557       }
4558
4559     return false;
4560   }
4561   //================================================================================
4562   /*!
4563    * \brief Tries to create a pentahedron
4564    */
4565   bool Hexahedron::addPenta()
4566   {
4567     // find a base triangular face
4568     int iTri = -1;
4569     for ( int iF = 0; iF < 5 && iTri < 0; ++iF )
4570       if ( _polygons[ iF ]._links.size() == 3 )
4571         iTri = iF;
4572     if ( iTri < 0 ) return false;
4573
4574     // find nodes
4575     _Node* nodes[6];
4576     int nbN = 0;
4577     for ( int iL = 0; iL < 3; ++iL )
4578     {
4579       // a base node
4580       nodes[iL] = _polygons[ iTri ]._links[iL].FirstNode();
4581       ++nbN;
4582
4583       // find a top node above the base node
4584       _Link* link = _polygons[ iTri ]._links[iL]._link;
4585       if ( !link->_faces[0] || !link->_faces[1] )
4586         return debugDumpLink( link );
4587       // a quadrangle sharing <link> with a base triangle
4588       _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[ iTri ] )];
4589       if ( quad->_links.size() != 4 ) return false;
4590       for ( int i = 0; i < 4; ++i )
4591         if ( quad->_links[i]._link == link )
4592         {
4593           // 1st node of a link opposite to <link> in <quad>
4594           nodes[iL+3] = quad->_links[(i+2)%4].FirstNode();
4595           ++nbN;
4596           break;
4597         }
4598     }
4599     if ( nbN == 6 )
4600       _volumeDefs.Set( &nodes[0], 6 );
4601
4602     return ( nbN == 6 );
4603   }
4604   //================================================================================
4605   /*!
4606    * \brief Tries to create a pyramid
4607    */
4608   bool Hexahedron::addPyra()
4609   {
4610     // find a base quadrangle
4611     int iQuad = -1;
4612     for ( int iF = 0; iF < 5 && iQuad < 0; ++iF )
4613       if ( _polygons[ iF ]._links.size() == 4 )
4614         iQuad = iF;
4615     if ( iQuad < 0 ) return false;
4616
4617     // find nodes
4618     _Node* nodes[5];
4619     nodes[0] = _polygons[iQuad]._links[0].FirstNode();
4620     nodes[1] = _polygons[iQuad]._links[1].FirstNode();
4621     nodes[2] = _polygons[iQuad]._links[2].FirstNode();
4622     nodes[3] = _polygons[iQuad]._links[3].FirstNode();
4623
4624     _Link* link = _polygons[iQuad]._links[0]._link;
4625     if ( !link->_faces[0] || !link->_faces[1] )
4626       return debugDumpLink( link );
4627
4628     // a triangle sharing <link> with a base quadrangle
4629     _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[ iQuad ] )];
4630     if ( tria->_links.size() != 3 ) return false;
4631     for ( int i = 0; i < 3; ++i )
4632       if ( tria->_links[i]._link == link )
4633       {
4634         nodes[4] = tria->_links[(i+1)%3].LastNode();
4635         _volumeDefs.Set( &nodes[0], 5 );
4636         return true;
4637       }
4638
4639     return false;
4640   }
4641   //================================================================================
4642   /*!
4643    * \brief Dump a link and return \c false
4644    */
4645   bool Hexahedron::debugDumpLink( Hexahedron::_Link* link )
4646   {
4647 #ifdef _DEBUG_
4648     gp_Pnt p1 = link->_nodes[0]->Point(), p2 = link->_nodes[1]->Point();
4649     cout << "BUG: not shared link. IKJ = ( "<< _i << " " << _j << " " << _k << " )" << endl
4650          << "n1 (" << p1.X() << ", "<< p1.Y() << ", "<< p1.Z() << " )" << endl
4651          << "n2 (" << p2.X() << ", "<< p2.Y() << ", "<< p2.Z() << " )" << endl;
4652 #endif
4653     return false;
4654   }
4655   //================================================================================
4656   /*!
4657    * \brief Classify a point by grid parameters
4658    */
4659   bool Hexahedron::isOutParam(const double uvw[3]) const
4660   {
4661     return (( _grid->_coords[0][ _i   ] - _grid->_tol > uvw[0] ) ||
4662             ( _grid->_coords[0][ _i+1 ] + _grid->_tol < uvw[0] ) ||
4663             ( _grid->_coords[1][ _j   ] - _grid->_tol > uvw[1] ) ||
4664             ( _grid->_coords[1][ _j+1 ] + _grid->_tol < uvw[1] ) ||
4665             ( _grid->_coords[2][ _k   ] - _grid->_tol > uvw[2] ) ||
4666             ( _grid->_coords[2][ _k+1 ] + _grid->_tol < uvw[2] ));
4667   }
4668   //================================================================================
4669   /*!
4670    * \brief Divide a polygon into triangles and modify accordingly an adjacent polyhedron
4671    */
4672   void splitPolygon( const SMDS_MeshElement*         polygon,
4673                      SMDS_VolumeTool &               volume,
4674                      const int                       facetIndex,
4675                      const TGeomID                   faceID,
4676                      const TGeomID                   solidID,
4677                      SMESH_MeshEditor::ElemFeatures& face,
4678                      SMESH_MeshEditor&               editor,
4679                      const bool                      reinitVolume)
4680   {
4681     SMESH_MeshAlgos::Triangulate divider(/*optimize=*/false);
4682     int nbTrias = divider.GetTriangles( polygon, face.myNodes );
4683     face.myNodes.resize( nbTrias * 3 );
4684
4685     SMESH_MeshEditor::ElemFeatures newVolumeDef;
4686     newVolumeDef.Init( volume.Element() );
4687     newVolumeDef.SetID( volume.Element()->GetID() );
4688
4689     newVolumeDef.myPolyhedQuantities.reserve( volume.NbFaces() + nbTrias );
4690     newVolumeDef.myNodes.reserve( volume.NbNodes() + nbTrias * 3 );
4691
4692     SMESHDS_Mesh* meshDS = editor.GetMeshDS();
4693     SMDS_MeshElement* newTriangle;
4694     for ( int iF = 0, nF = volume.NbFaces(); iF < nF; iF++ )
4695     {
4696       if ( iF == facetIndex )
4697       {
4698         newVolumeDef.myPolyhedQuantities.push_back( 3 );
4699         newVolumeDef.myNodes.insert( newVolumeDef.myNodes.end(),
4700                                      face.myNodes.begin(),
4701                                      face.myNodes.begin() + 3 );
4702         meshDS->RemoveFreeElement( polygon, 0, false );
4703         newTriangle = meshDS->AddFace( face.myNodes[0], face.myNodes[1], face.myNodes[2] );
4704         meshDS->SetMeshElementOnShape( newTriangle, faceID );
4705       }
4706       else
4707       {
4708         const SMDS_MeshNode** nn = volume.GetFaceNodes( iF );
4709         const size_t nbFaceNodes = volume.NbFaceNodes ( iF );
4710         newVolumeDef.myPolyhedQuantities.push_back( nbFaceNodes );
4711         newVolumeDef.myNodes.insert( newVolumeDef.myNodes.end(), nn, nn + nbFaceNodes );
4712       }
4713     }
4714
4715     for ( size_t iN = 3; iN < face.myNodes.size(); iN += 3 )
4716     {
4717       newVolumeDef.myPolyhedQuantities.push_back( 3 );
4718       newVolumeDef.myNodes.insert( newVolumeDef.myNodes.end(),
4719                                    face.myNodes.begin() + iN,
4720                                    face.myNodes.begin() + iN + 3 );
4721       newTriangle = meshDS->AddFace( face.myNodes[iN], face.myNodes[iN+1], face.myNodes[iN+2] );
4722       meshDS->SetMeshElementOnShape( newTriangle, faceID );
4723     }
4724
4725     meshDS->RemoveFreeElement( volume.Element(), 0, false );
4726     SMDS_MeshElement* newVolume = editor.AddElement( newVolumeDef.myNodes, newVolumeDef );
4727     meshDS->SetMeshElementOnShape( newVolume, solidID );
4728
4729     if ( reinitVolume )
4730     {
4731       volume.Set( 0 );
4732       volume.Set( newVolume );
4733     }
4734     return;
4735   }
4736   //================================================================================
4737   /*!
4738    * \brief Create mesh faces at free facets
4739    */
4740   void Hexahedron::addFaces( SMESH_MesherHelper&                       helper,
4741                              const vector< const SMDS_MeshElement* > & boundaryVolumes )
4742   {
4743     if ( !_grid->_toCreateFaces )
4744       return;
4745
4746     SMDS_VolumeTool vTool;
4747     vector<int> bndFacets;
4748     SMESH_MeshEditor editor( helper.GetMesh() );
4749     SMESH_MeshEditor::ElemFeatures face( SMDSAbs_Face );
4750     SMESHDS_Mesh* meshDS = helper.GetMeshDS();
4751
4752     // check if there are internal or shared FACEs
4753     bool hasInternal = ( !_grid->_geometry.IsOneSolid() ||
4754                          _grid->_geometry._soleSolid.HasInternalFaces() );
4755
4756     for ( size_t iV = 0; iV < boundaryVolumes.size(); ++iV )
4757     {
4758       if ( !vTool.Set( boundaryVolumes[ iV ]))
4759         continue;
4760
4761       TGeomID solidID = vTool.Element()->GetShapeID();
4762       Solid *   solid = _grid->GetOneOfSolids( solidID );
4763
4764       // find boundary facets
4765
4766       bndFacets.clear();
4767       for ( int iF = 0, n = vTool.NbFaces(); iF < n; iF++ )
4768       {
4769         bool isBoundary = vTool.IsFreeFace( iF );
4770         if ( isBoundary )
4771         {
4772           bndFacets.push_back( iF );
4773         }
4774         else if ( hasInternal )
4775         {
4776           // check if all nodes are on internal/shared FACEs
4777           isBoundary = true;
4778           const SMDS_MeshNode** nn = vTool.GetFaceNodes( iF );
4779           const size_t nbFaceNodes = vTool.NbFaceNodes ( iF );
4780           for ( size_t iN = 0; iN < nbFaceNodes &&  isBoundary; ++iN )
4781             isBoundary = ( nn[ iN ]->GetShapeID() != solidID );
4782           if ( isBoundary )
4783             bndFacets.push_back( -( iF+1 )); // !!! minus ==> to check the FACE
4784         }
4785       }
4786       if ( bndFacets.empty() )
4787         continue;
4788
4789       // create faces
4790
4791       if ( !vTool.IsPoly() )
4792         vTool.SetExternalNormal();
4793       for ( size_t i = 0; i < bndFacets.size(); ++i ) // loop on boundary facets
4794       {
4795         const bool    isBoundary = ( bndFacets[i] >= 0 );
4796         const int         iFacet = isBoundary ? bndFacets[i] : -bndFacets[i]-1;
4797         const SMDS_MeshNode** nn = vTool.GetFaceNodes( iFacet );
4798         const size_t nbFaceNodes = vTool.NbFaceNodes ( iFacet );
4799         face.myNodes.assign( nn, nn + nbFaceNodes );
4800
4801         TGeomID faceID = 0;
4802         const SMDS_MeshElement* existFace = 0, *newFace = 0;
4803
4804         if (( existFace = meshDS->FindElement( face.myNodes, SMDSAbs_Face )))
4805         {
4806           if ( existFace->isMarked() )
4807             continue; // created by this method
4808           faceID = existFace->GetShapeID();
4809         }
4810         else
4811         {
4812           // look for a supporting FACE
4813           for ( size_t iN = 0; iN < nbFaceNodes &&  !faceID; ++iN ) // look for a node on FACE
4814           {
4815             if ( nn[ iN ]->GetPosition()->GetDim() == 2 )
4816               faceID = nn[ iN ]->GetShapeID();
4817           }
4818           for ( size_t iN = 0; iN < nbFaceNodes &&  !faceID; ++iN )
4819           {
4820             // look for a father FACE of EDGEs and VERTEXes
4821             const TopoDS_Shape& s1 = _grid->Shape( nn[ iN   ]->GetShapeID() );
4822             const TopoDS_Shape& s2 = _grid->Shape( nn[ iN+1 ]->GetShapeID() );
4823             if ( s1 != s2 && s1.ShapeType() == TopAbs_EDGE && s2.ShapeType() == TopAbs_EDGE )
4824             {
4825               TopoDS_Shape f = helper.GetCommonAncestor( s1, s2, *helper.GetMesh(), TopAbs_FACE );
4826               if ( !f.IsNull() )
4827                 faceID = _grid->ShapeID( f );
4828             }
4829           }
4830
4831           bool toCheckFace = faceID && (( !isBoundary ) ||
4832                                         ( hasInternal && _grid->_toUseThresholdForInternalFaces ));
4833           if ( toCheckFace ) // check if all nodes are on the found FACE
4834           {
4835             SMESH_subMesh* faceSM = helper.GetMesh()->GetSubMeshContaining( faceID );
4836             for ( size_t iN = 0; iN < nbFaceNodes &&  faceID; ++iN )
4837             {
4838               TGeomID subID = nn[ iN ]->GetShapeID();
4839               if ( subID != faceID && !faceSM->DependsOn( subID ))
4840                 faceID = 0;
4841             }
4842             if ( !faceID && !isBoundary )
4843               continue;
4844           }
4845         }
4846         // orient a new face according to supporting FACE orientation in shape_to_mesh
4847         if ( !solid->IsOutsideOriented( faceID ))
4848         {
4849           if ( existFace )
4850             editor.Reorient( existFace );
4851           else
4852             std::reverse( face.myNodes.begin(), face.myNodes.end() );
4853         }
4854
4855         if ( ! ( newFace = existFace ))
4856         {
4857           face.SetPoly( nbFaceNodes > 4 );
4858           newFace = editor.AddElement( face.myNodes, face );
4859           if ( !newFace )
4860             continue;
4861           newFace->setIsMarked( true ); // to distinguish from face created in getBoundaryElems()
4862         }
4863
4864         if ( faceID && _grid->IsBoundaryFace( faceID )) // face is not shared
4865         {
4866           // set newFace to the found FACE provided that it fully lies on the FACE
4867           for ( size_t iN = 0; iN < nbFaceNodes &&  faceID; ++iN )
4868             if ( nn[iN]->GetShapeID() == solidID )
4869             {
4870               if ( existFace )
4871                 meshDS->UnSetMeshElementOnShape( existFace, _grid->Shape( faceID ));
4872               faceID = 0;
4873             }
4874         }
4875
4876         // split a polygon that will be used by other 3D algorithm
4877         if ( faceID && nbFaceNodes > 4 &&
4878              !_grid->IsInternal( faceID ) &&
4879              !_grid->IsShared( faceID ) &&
4880              !_grid->IsBoundaryFace( faceID ))
4881         {
4882           splitPolygon( newFace, vTool, iFacet, faceID, solidID,
4883                         face, editor, i+1 < bndFacets.size() );
4884         }
4885         else
4886         {
4887           if ( faceID )
4888             meshDS->SetMeshElementOnShape( newFace, faceID );
4889           else
4890             meshDS->SetMeshElementOnShape( newFace, solidID );
4891         }
4892       } // loop on bndFacets
4893     } // loop on boundaryVolumes
4894
4895
4896     // Orient coherently mesh faces on INTERNAL FACEs
4897
4898     if ( hasInternal )
4899     {
4900       TopExp_Explorer faceExp( _grid->_geometry._mainShape, TopAbs_FACE );
4901       for ( ; faceExp.More(); faceExp.Next() )
4902       {
4903         if ( faceExp.Current().Orientation() != TopAbs_INTERNAL )
4904           continue;
4905
4906         SMESHDS_SubMesh* sm = meshDS->MeshElements( faceExp.Current() );
4907         if ( !sm ) continue;
4908
4909         TIDSortedElemSet facesToOrient;
4910         for ( SMDS_ElemIteratorPtr fIt = sm->GetElements(); fIt->more(); )
4911           facesToOrient.insert( facesToOrient.end(), fIt->next() );
4912         if ( facesToOrient.size() < 2 )
4913           continue;
4914
4915         gp_Dir direction(1,0,0);
4916         const SMDS_MeshElement* anyFace = *facesToOrient.begin();
4917         editor.Reorient2D( facesToOrient, direction, anyFace );
4918       }
4919     }
4920     return;
4921   }
4922
4923   //================================================================================
4924   /*!
4925    * \brief Create mesh segments.
4926    */
4927   void Hexahedron::addSegments( SMESH_MesherHelper&                      helper,
4928                                 const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap )
4929   {
4930     SMESHDS_Mesh* mesh = helper.GetMeshDS();
4931
4932     std::vector<const SMDS_MeshNode*> nodes;
4933     std::vector<const SMDS_MeshElement *> elems;
4934     map< TGeomID, vector< TGeomID > >::const_iterator e2ff = edge2faceIDsMap.begin();
4935     for ( ; e2ff != edge2faceIDsMap.end(); ++e2ff )
4936     {
4937       const TopoDS_Edge& edge = TopoDS::Edge( _grid->Shape( e2ff->first ));
4938       const TopoDS_Face& face = TopoDS::Face( _grid->Shape( e2ff->second[0] ));
4939       StdMeshers_FaceSide side( face, edge, helper.GetMesh(), /*isFwd=*/true, /*skipMed=*/true );
4940       nodes = side.GetOrderedNodes();
4941
4942       elems.clear();
4943       if ( nodes.size() == 2 )
4944         // check that there is an element connecting two nodes
4945         if ( !mesh->GetElementsByNodes( nodes, elems ))
4946           continue;
4947
4948       for ( size_t i = 1; i < nodes.size(); i++ )
4949       {
4950         SMDS_MeshElement* segment = mesh->AddEdge( nodes[i-1], nodes[i] );
4951         mesh->SetMeshElementOnShape( segment, e2ff->first );
4952       }
4953     }
4954     return;
4955   }
4956
4957   //================================================================================
4958   /*!
4959    * \brief Return created volumes and volumes that can have free facet because of
4960    *        skipped small volume. Also create mesh faces on free facets
4961    *        of adjacent not-cut volumes id the result volume is too small.
4962    */
4963   void Hexahedron::getBoundaryElems( vector< const SMDS_MeshElement* > & boundaryElems )
4964   {
4965     if ( _hasTooSmall /*|| _volumeDefs.IsEmpty()*/ )
4966     {
4967       // create faces around a missing small volume
4968       TGeomID faceID = 0;
4969       SMESH_MeshEditor editor( _grid->_helper->GetMesh() );
4970       SMESH_MeshEditor::ElemFeatures polygon( SMDSAbs_Face );
4971       SMESHDS_Mesh* meshDS = _grid->_helper->GetMeshDS();
4972       std::vector<const SMDS_MeshElement *> adjVolumes(2);
4973       for ( size_t iF = 0; iF < _polygons.size(); ++iF )
4974       {
4975         const size_t nbLinks = _polygons[ iF ]._links.size();
4976         if ( nbLinks != 4 ) continue;
4977         polygon.myNodes.resize( nbLinks );
4978         polygon.myNodes.back() = 0;
4979         for ( size_t iL = 0, iN = nbLinks - 1; iL < nbLinks; ++iL, --iN )
4980           if ( ! ( polygon.myNodes[iN] = _polygons[ iF ]._links[ iL ].FirstNode()->Node() ))
4981             break;
4982         if ( !polygon.myNodes.back() )
4983           continue;
4984
4985         meshDS->GetElementsByNodes( polygon.myNodes, adjVolumes, SMDSAbs_Volume );
4986         if ( adjVolumes.size() != 1 )
4987           continue;
4988         if ( !adjVolumes[0]->isMarked() )
4989         {
4990           boundaryElems.push_back( adjVolumes[0] );
4991           adjVolumes[0]->setIsMarked( true );
4992         }
4993
4994         bool sameShape = true;
4995         TGeomID shapeID = polygon.myNodes[0]->GetShapeID();
4996         for ( size_t i = 1; i < polygon.myNodes.size() && sameShape; ++i )
4997           sameShape = ( shapeID == polygon.myNodes[i]->GetShapeID() );
4998
4999         if ( !sameShape || !_grid->IsSolid( shapeID ))
5000           continue; // some of shapes must be FACE
5001
5002         if ( !faceID )
5003         {
5004           faceID = getAnyFace();
5005           if ( !faceID )
5006             break;
5007           if ( _grid->IsInternal( faceID ) ||
5008                _grid->IsShared( faceID ) ||
5009                _grid->IsBoundaryFace( faceID ))
5010             break; // create only if a new face will be used by other 3D algo
5011         }
5012
5013         Solid * solid = _grid->GetOneOfSolids( adjVolumes[0]->GetShapeID() );
5014         if ( !solid->IsOutsideOriented( faceID ))
5015           std::reverse( polygon.myNodes.begin(), polygon.myNodes.end() );
5016
5017         //polygon.SetPoly( polygon.myNodes.size() > 4 );
5018         const SMDS_MeshElement* newFace = editor.AddElement( polygon.myNodes, polygon );
5019         meshDS->SetMeshElementOnShape( newFace, faceID );
5020       }
5021     }
5022
5023     // return created volumes
5024     for ( _volumeDef* volDef = &_volumeDefs; volDef; volDef = volDef->_next )
5025     {
5026       if ( volDef->_volume && !volDef->_volume->isMarked() )
5027       {
5028         volDef->_volume->setIsMarked( true );
5029         boundaryElems.push_back( volDef->_volume );
5030
5031         if ( _grid->IsToCheckNodePos() ) // un-mark nodes marked in addVolumes()
5032           for ( size_t iN = 0; iN < volDef->_nodes.size(); ++iN )
5033             volDef->_nodes[iN].Node()->setIsMarked( false );
5034       }
5035     }
5036   }
5037
5038   //================================================================================
5039   /*!
5040    * \brief Set to _hexLinks a next portion of splits located on one side of INTERNAL FACEs
5041    */
5042   bool Hexahedron::_SplitIterator::Next()
5043   {
5044     if ( _iterationNb > 0 )
5045       // count used splits
5046       for ( size_t i = 0; i < _splits.size(); ++i )
5047       {
5048         if ( _splits[i]._iCheckIteration == _iterationNb )
5049         {
5050           _splits[i]._isUsed = _splits[i]._checkedSplit->_faces[1];
5051           _nbUsed += _splits[i]._isUsed;
5052         }
5053         if ( !More() )
5054           return false;
5055       }
5056
5057     ++_iterationNb;
5058
5059     bool toTestUsed = ( _nbChecked >= _splits.size() );
5060     if ( toTestUsed )
5061     {
5062       // all splits are checked; find all not used splits
5063       for ( size_t i = 0; i < _splits.size(); ++i )
5064         if ( !_splits[i].IsCheckedOrUsed( toTestUsed ))
5065           _splits[i]._iCheckIteration = _iterationNb;
5066
5067       _nbUsed = _splits.size(); // to stop iteration
5068     }
5069     else
5070     {
5071       // get any not used/checked split to start from
5072       _freeNodes.clear();
5073       for ( size_t i = 0; i < _splits.size(); ++i )
5074       {
5075         if ( !_splits[i].IsCheckedOrUsed( toTestUsed ))
5076         {
5077           _freeNodes.push_back( _splits[i]._nodes[0] );
5078           _freeNodes.push_back( _splits[i]._nodes[1] );
5079           _splits[i]._iCheckIteration = _iterationNb;
5080           break;
5081         }
5082       }
5083       // find splits connected to the start one via _freeNodes
5084       for ( size_t iN = 0; iN < _freeNodes.size(); ++iN )
5085       {
5086         for ( size_t iS = 0; iS < _splits.size(); ++iS )
5087         {
5088           if ( _splits[iS].IsCheckedOrUsed( toTestUsed ))
5089             continue;
5090           int iN2 = -1;
5091           if (      _freeNodes[iN] == _splits[iS]._nodes[0] )
5092             iN2 = 1;
5093           else if ( _freeNodes[iN] == _splits[iS]._nodes[1] )
5094             iN2 = 0;
5095           else
5096             continue;
5097           if ( _freeNodes[iN]->_isInternalFlags > 0 )
5098           {
5099             if ( _splits[iS]._nodes[ iN2 ]->_isInternalFlags == 0 )
5100               continue;
5101             if ( !_splits[iS]._nodes[ iN2 ]->IsLinked( _freeNodes[iN]->_intPoint ))
5102               continue;
5103           }
5104           _splits[iS]._iCheckIteration = _iterationNb;
5105           _freeNodes.push_back( _splits[iS]._nodes[ iN2 ]);
5106         }
5107       }
5108     }
5109     // set splits to hex links
5110
5111     for ( int iL = 0; iL < 12; ++iL )
5112       _hexLinks[ iL ]._splits.clear();
5113
5114     _Link split;
5115     for ( size_t i = 0; i < _splits.size(); ++i )
5116     {
5117       if ( _splits[i]._iCheckIteration == _iterationNb )
5118       {
5119         split._nodes[0] = _splits[i]._nodes[0];
5120         split._nodes[1] = _splits[i]._nodes[1];
5121         _Link & hexLink = _hexLinks[ _splits[i]._linkID ];
5122         hexLink._splits.push_back( split );
5123         _splits[i]._checkedSplit = & hexLink._splits.back();
5124         ++_nbChecked;
5125       }
5126     }
5127     return More();
5128   }
5129
5130   //================================================================================
5131   /*!
5132    * \brief computes exact bounding box with axes parallel to given ones
5133    */
5134   //================================================================================
5135
5136   void getExactBndBox( const vector< TopoDS_Shape >& faceVec,
5137                        const double*                 axesDirs,
5138                        Bnd_Box&                      shapeBox )
5139   {
5140     BRep_Builder b;
5141     TopoDS_Compound allFacesComp;
5142     b.MakeCompound( allFacesComp );
5143     for ( size_t iF = 0; iF < faceVec.size(); ++iF )
5144       b.Add( allFacesComp, faceVec[ iF ] );
5145
5146     double sP[6]; // aXmin, aYmin, aZmin, aXmax, aYmax, aZmax
5147     shapeBox.Get(sP[0],sP[1],sP[2],sP[3],sP[4],sP[5]);
5148     double farDist = 0;
5149     for ( int i = 0; i < 6; ++i )
5150       farDist = Max( farDist, 10 * sP[i] );
5151
5152     gp_XYZ axis[3] = { gp_XYZ( axesDirs[0], axesDirs[1], axesDirs[2] ),
5153                        gp_XYZ( axesDirs[3], axesDirs[4], axesDirs[5] ),
5154                        gp_XYZ( axesDirs[6], axesDirs[7], axesDirs[8] ) };
5155     axis[0].Normalize();
5156     axis[1].Normalize();
5157     axis[2].Normalize();
5158
5159     gp_Mat basis( axis[0], axis[1], axis[2] );
5160     gp_Mat bi = basis.Inverted();
5161
5162     gp_Pnt pMin, pMax;
5163     for ( int iDir = 0; iDir < 3; ++iDir )
5164     {
5165       gp_XYZ axis0 = axis[ iDir ];
5166       gp_XYZ axis1 = axis[ ( iDir + 1 ) % 3 ];
5167       gp_XYZ axis2 = axis[ ( iDir + 2 ) % 3 ];
5168       for ( int isMax = 0; isMax < 2; ++isMax )
5169       {
5170         double shift = isMax ? farDist : -farDist;
5171         gp_XYZ orig = shift * axis0;
5172         gp_XYZ norm = axis1 ^ axis2;
5173         gp_Pln pln( orig, norm );
5174         norm = pln.Axis().Direction().XYZ();
5175         BRepBuilderAPI_MakeFace plane( pln, -farDist, farDist, -farDist, farDist );
5176
5177         gp_Pnt& pAxis = isMax ? pMax : pMin;
5178         gp_Pnt pPlane, pFaces;
5179         double dist = GEOMUtils::GetMinDistance( plane, allFacesComp, pPlane, pFaces );
5180         if ( dist < 0 )
5181         {
5182           Bnd_B3d bb;
5183           gp_XYZ corner;
5184           for ( int i = 0; i < 2; ++i ) {
5185             corner.SetCoord( 1, sP[ i*3 ]);
5186             for ( int j = 0; j < 2; ++j ) {
5187               corner.SetCoord( 2, sP[ i*3 + 1 ]);
5188               for ( int k = 0; k < 2; ++k )
5189               {
5190                 corner.SetCoord( 3, sP[ i*3 + 2 ]);
5191                 corner *= bi;
5192                 bb.Add( corner );
5193               }
5194             }
5195           }
5196           corner = isMax ? bb.CornerMax() : bb.CornerMin();
5197           pAxis.SetCoord( iDir+1, corner.Coord( iDir+1 ));
5198         }
5199         else
5200         {
5201           gp_XYZ pf = pFaces.XYZ() * bi;
5202           pAxis.SetCoord( iDir+1, pf.Coord( iDir+1 ) );
5203         }
5204       }
5205     } // loop on 3 axes
5206
5207     shapeBox.SetVoid();
5208     shapeBox.Add( pMin );
5209     shapeBox.Add( pMax );
5210
5211     return;
5212   }
5213
5214 } // namespace
5215
5216 //=============================================================================
5217 /*!
5218  * \brief Generates 3D structured Cartesian mesh in the internal part of
5219  * solid shapes and polyhedral volumes near the shape boundary.
5220  *  \param theMesh - mesh to fill in
5221  *  \param theShape - a compound of all SOLIDs to mesh
5222  *  \retval bool - true in case of success
5223  */
5224 //=============================================================================
5225
5226 bool StdMeshers_Cartesian_3D::Compute(SMESH_Mesh &         theMesh,
5227                                       const TopoDS_Shape & theShape)
5228 {
5229   // The algorithm generates the mesh in following steps:
5230
5231   // 1) Intersection of grid lines with the geometry boundary.
5232   // This step allows to find out if a given node of the initial grid is
5233   // inside or outside the geometry.
5234
5235   // 2) For each cell of the grid, check how many of it's nodes are outside
5236   // of the geometry boundary. Depending on a result of this check
5237   // - skip a cell, if all it's nodes are outside
5238   // - skip a cell, if it is too small according to the size threshold
5239   // - add a hexahedron in the mesh, if all nodes are inside
5240   // - add a polyhedron in the mesh, if some nodes are inside and some outside
5241
5242   _computeCanceled = false;
5243
5244   SMESH_MesherHelper helper( theMesh );
5245   SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
5246
5247   try
5248   {
5249     Grid grid;
5250     grid._helper                         = &helper;
5251     grid._toAddEdges                     = _hyp->GetToAddEdges();
5252     grid._toCreateFaces                  = _hyp->GetToCreateFaces();
5253     grid._toConsiderInternalFaces        = _hyp->GetToConsiderInternalFaces();
5254     grid._toUseThresholdForInternalFaces = _hyp->GetToUseThresholdForInternalFaces();
5255     grid._sizeThreshold                  = _hyp->GetSizeThreshold();
5256     grid.InitGeometry( theShape );
5257
5258     vector< TopoDS_Shape > faceVec;
5259     {
5260       TopTools_MapOfShape faceMap;
5261       TopExp_Explorer fExp;
5262       for ( fExp.Init( theShape, TopAbs_FACE ); fExp.More(); fExp.Next() )
5263       {
5264         bool isNewFace = faceMap.Add( fExp.Current() );
5265         if ( !grid._toConsiderInternalFaces )
5266           if ( !isNewFace || fExp.Current().Orientation() == TopAbs_INTERNAL )
5267             // remove an internal face
5268             faceMap.Remove( fExp.Current() );
5269       }
5270       faceVec.reserve( faceMap.Extent() );
5271       faceVec.assign( faceMap.cbegin(), faceMap.cend() );
5272     }
5273     vector<FaceGridIntersector> facesItersectors( faceVec.size() );
5274     Bnd_Box shapeBox;
5275     for ( size_t i = 0; i < faceVec.size(); ++i )
5276     {
5277       facesItersectors[i]._face   = TopoDS::Face( faceVec[i] );
5278       facesItersectors[i]._faceID = grid.ShapeID( faceVec[i] );
5279       facesItersectors[i]._grid   = &grid;
5280       shapeBox.Add( facesItersectors[i].GetFaceBndBox() );
5281     }
5282     getExactBndBox( faceVec, _hyp->GetAxisDirs(), shapeBox );
5283
5284
5285     vector<double> xCoords, yCoords, zCoords;
5286     _hyp->GetCoordinates( xCoords, yCoords, zCoords, shapeBox );
5287
5288     grid.SetCoordinates( xCoords, yCoords, zCoords, _hyp->GetAxisDirs(), shapeBox );
5289
5290     if ( _computeCanceled ) return false;
5291
5292 #ifdef WITH_TBB
5293     { // copy partner faces and curves of not thread-safe types
5294       set< const Standard_Transient* > tshapes;
5295       BRepBuilderAPI_Copy copier;
5296       for ( size_t i = 0; i < facesItersectors.size(); ++i )
5297       {
5298         if ( !facesItersectors[i].IsThreadSafe( tshapes ))
5299         {
5300           copier.Perform( facesItersectors[i]._face );
5301           facesItersectors[i]._face = TopoDS::Face( copier );
5302         }
5303       }
5304     }
5305     // Intersection of grid lines with the geometry boundary.
5306     tbb::parallel_for ( tbb::blocked_range<size_t>( 0, facesItersectors.size() ),
5307                         ParallelIntersector( facesItersectors ),
5308                         tbb::simple_partitioner());
5309 #else
5310     for ( size_t i = 0; i < facesItersectors.size(); ++i )
5311       facesItersectors[i].Intersect();
5312 #endif
5313
5314     // put intersection points onto the GridLine's; this is done after intersection
5315     // to avoid contention of facesItersectors for writing into the same GridLine
5316     // in case of parallel work of facesItersectors
5317     for ( size_t i = 0; i < facesItersectors.size(); ++i )
5318       facesItersectors[i].StoreIntersections();
5319
5320     if ( _computeCanceled ) return false;
5321
5322     // create nodes on the geometry
5323     grid.ComputeNodes( helper );
5324
5325     if ( _computeCanceled ) return false;
5326
5327     // get EDGEs to take into account
5328     map< TGeomID, vector< TGeomID > > edge2faceIDsMap;
5329     grid.GetEdgesToImplement( edge2faceIDsMap, theShape, faceVec );
5330
5331     // create volume elements
5332     Hexahedron hex( &grid );
5333     int nbAdded = hex.MakeElements( helper, edge2faceIDsMap );
5334
5335     if ( nbAdded > 0 )
5336     {
5337       if ( !grid._toConsiderInternalFaces )
5338       {
5339         // make all SOLIDs computed
5340         TopExp_Explorer solidExp( theShape, TopAbs_SOLID );
5341         if ( SMESHDS_SubMesh* sm1 = meshDS->MeshElements( solidExp.Current()) )
5342         {
5343           SMDS_ElemIteratorPtr volIt = sm1->GetElements();
5344           for ( ; solidExp.More() && volIt->more(); solidExp.Next() )
5345           {
5346             const SMDS_MeshElement* vol = volIt->next();
5347             sm1->RemoveElement( vol );
5348             meshDS->SetMeshElementOnShape( vol, solidExp.Current() );
5349           }
5350         }
5351       }
5352       // make other sub-shapes computed
5353       setSubmeshesComputed( theMesh, theShape );
5354     }
5355
5356     // remove free nodes
5357     //if ( SMESHDS_SubMesh * smDS = meshDS->MeshElements( helper.GetSubShapeID() ))
5358     {
5359       std::vector< const SMDS_MeshNode* > nodesToRemove;
5360       // get intersection nodes
5361       for ( int iDir = 0; iDir < 3; ++iDir )
5362       {
5363         vector< GridLine >& lines = grid._lines[ iDir ];
5364         for ( size_t i = 0; i < lines.size(); ++i )
5365         {
5366           multiset< F_IntersectPoint >::iterator ip = lines[i]._intPoints.begin();
5367           for ( ; ip != lines[i]._intPoints.end(); ++ip )
5368             if ( ip->_node && ip->_node->NbInverseElements() == 0 && !ip->_node->isMarked() )
5369             {
5370               nodesToRemove.push_back( ip->_node );
5371               ip->_node->setIsMarked( true );
5372             }
5373         }
5374       }
5375       // get grid nodes
5376       for ( size_t i = 0; i < grid._nodes.size(); ++i )
5377         if ( grid._nodes[i] && grid._nodes[i]->NbInverseElements() == 0 &&
5378              !grid._nodes[i]->isMarked() )
5379         {
5380           nodesToRemove.push_back( grid._nodes[i] );
5381           grid._nodes[i]->setIsMarked( true );
5382         }
5383
5384       // do remove
5385       for ( size_t i = 0; i < nodesToRemove.size(); ++i )
5386         meshDS->RemoveFreeNode( nodesToRemove[i], /*smD=*/0, /*fromGroups=*/false );
5387     }
5388
5389     return nbAdded;
5390
5391   }
5392   // SMESH_ComputeError is not caught at SMESH_submesh level for an unknown reason
5393   catch ( SMESH_ComputeError& e)
5394   {
5395     return error( SMESH_ComputeErrorPtr( new SMESH_ComputeError( e )));
5396   }
5397   return false;
5398 }
5399
5400 //=============================================================================
5401 /*!
5402  *  Evaluate
5403  */
5404 //=============================================================================
5405
5406 bool StdMeshers_Cartesian_3D::Evaluate(SMESH_Mesh &         theMesh,
5407                                        const TopoDS_Shape & theShape,
5408                                        MapShapeNbElems&     theResMap)
5409 {
5410   // TODO
5411 //   std::vector<int> aResVec(SMDSEntity_Last);
5412 //   for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
5413 //   if(IsQuadratic) {
5414 //     aResVec[SMDSEntity_Quad_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
5415 //     int nb1d_face0_int = ( nb2d_face0*4 - nb1d ) / 2;
5416 //     aResVec[SMDSEntity_Node] = nb0d_face0 * ( 2*nb2d/nb1d - 1 ) - nb1d_face0_int * nb2d/nb1d;
5417 //   }
5418 //   else {
5419 //     aResVec[SMDSEntity_Node] = nb0d_face0 * ( nb2d/nb1d - 1 );
5420 //     aResVec[SMDSEntity_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
5421 //   }
5422 //   SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
5423 //   aResMap.insert(std::make_pair(sm,aResVec));
5424
5425   return true;
5426 }
5427
5428 //=============================================================================
5429 namespace
5430 {
5431   /*!
5432    * \brief Event listener setting/unsetting _alwaysComputed flag to
5433    *        submeshes of inferior levels to prevent their computing
5434    */
5435   struct _EventListener : public SMESH_subMeshEventListener
5436   {
5437     string _algoName;
5438
5439     _EventListener(const string& algoName):
5440       SMESH_subMeshEventListener(/*isDeletable=*/true,"StdMeshers_Cartesian_3D::_EventListener"),
5441       _algoName(algoName)
5442     {}
5443     // --------------------------------------------------------------------------------
5444     // setting/unsetting _alwaysComputed flag to submeshes of inferior levels
5445     //
5446     static void setAlwaysComputed( const bool     isComputed,
5447                                    SMESH_subMesh* subMeshOfSolid)
5448     {
5449       SMESH_subMeshIteratorPtr smIt =
5450         subMeshOfSolid->getDependsOnIterator(/*includeSelf=*/false, /*complexShapeFirst=*/false);
5451       while ( smIt->more() )
5452       {
5453         SMESH_subMesh* sm = smIt->next();
5454         sm->SetIsAlwaysComputed( isComputed );
5455       }
5456       subMeshOfSolid->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
5457     }
5458
5459     // --------------------------------------------------------------------------------
5460     // unsetting _alwaysComputed flag if "Cartesian_3D" was removed
5461     //
5462     virtual void ProcessEvent(const int          event,
5463                               const int          eventType,
5464                               SMESH_subMesh*     subMeshOfSolid,
5465                               SMESH_subMeshEventListenerData* data,
5466                               const SMESH_Hypothesis*         hyp = 0)
5467     {
5468       if ( eventType == SMESH_subMesh::COMPUTE_EVENT )
5469       {
5470         setAlwaysComputed( subMeshOfSolid->GetComputeState() == SMESH_subMesh::COMPUTE_OK,
5471                            subMeshOfSolid );
5472       }
5473       else
5474       {
5475         SMESH_Algo* algo3D = subMeshOfSolid->GetAlgo();
5476         if ( !algo3D || _algoName != algo3D->GetName() )
5477           setAlwaysComputed( false, subMeshOfSolid );
5478       }
5479     }
5480
5481     // --------------------------------------------------------------------------------
5482     // set the event listener
5483     //
5484     static void SetOn( SMESH_subMesh* subMeshOfSolid, const string& algoName )
5485     {
5486       subMeshOfSolid->SetEventListener( new _EventListener( algoName ),
5487                                         /*data=*/0,
5488                                         subMeshOfSolid );
5489     }
5490
5491   }; // struct _EventListener
5492
5493 } // namespace
5494
5495 //================================================================================
5496 /*!
5497  * \brief Sets event listener to submeshes if necessary
5498  *  \param subMesh - submesh where algo is set
5499  * This method is called when a submesh gets HYP_OK algo_state.
5500  * After being set, event listener is notified on each event of a submesh.
5501  */
5502 //================================================================================
5503
5504 void StdMeshers_Cartesian_3D::SetEventListener(SMESH_subMesh* subMesh)
5505 {
5506   _EventListener::SetOn( subMesh, GetName() );
5507 }
5508
5509 //================================================================================
5510 /*!
5511  * \brief Set _alwaysComputed flag to submeshes of inferior levels to avoid their computing
5512  */
5513 //================================================================================
5514
5515 void StdMeshers_Cartesian_3D::setSubmeshesComputed(SMESH_Mesh&         theMesh,
5516                                                    const TopoDS_Shape& theShape)
5517 {
5518   for ( TopExp_Explorer soExp( theShape, TopAbs_SOLID ); soExp.More(); soExp.Next() )
5519     _EventListener::setAlwaysComputed( true, theMesh.GetSubMesh( soExp.Current() ));
5520 }