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