]> SALOME platform Git repositories - modules/smesh.git/blob - src/StdMeshers/StdMeshers_Cartesian_3D.cxx
Salome HOME
c392deeeaa41d3cc75f64671002a8a94f0b1842a
[modules/smesh.git] / src / StdMeshers / StdMeshers_Cartesian_3D.cxx
1 // Copyright (C) 2007-2024  CEA, EDF, 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 #include "StdMeshers_Cartesian_VL.hxx"
28 #include "StdMeshers_FaceSide.hxx"
29 #include "StdMeshers_ViscousLayers.hxx"
30
31 #include <ObjectPool.hxx>
32 #include <SMDS_LinearEdge.hxx>
33 #include <SMDS_MeshNode.hxx>
34 #include <SMDS_VolumeOfNodes.hxx>
35 #include <SMDS_VolumeTool.hxx>
36 #include <SMESHDS_Mesh.hxx>
37 #include <SMESH_Block.hxx>
38 #include <SMESH_Comment.hxx>
39 #include <SMESH_ControlsDef.hxx>
40 #include <SMESH_Mesh.hxx>
41 #include <SMESH_MeshAlgos.hxx>
42 #include <SMESH_MeshEditor.hxx>
43 #include <SMESH_MesherHelper.hxx>
44 #include <SMESH_subMesh.hxx>
45 #include <SMESH_subMeshEventListener.hxx>
46
47 #include <utilities.h>
48 #include <Utils_ExceptHandlers.hxx>
49
50 #include <GEOMUtils.hxx>
51
52 #include <BRepAdaptor_Curve.hxx>
53 #include <BRepAdaptor_Surface.hxx>
54 #include <BRepBndLib.hxx>
55 #include <BRepBuilderAPI_Copy.hxx>
56 #include <BRepBuilderAPI_MakeFace.hxx>
57 #include <BRepTools.hxx>
58 #include <BRepTopAdaptor_FClass2d.hxx>
59 #include <BRep_Builder.hxx>
60 #include <BRep_Tool.hxx>
61 #include <Bnd_B3d.hxx>
62 #include <Bnd_Box.hxx>
63 #include <ElSLib.hxx>
64 #include <GCPnts_UniformDeflection.hxx>
65 #include <Geom2d_BSplineCurve.hxx>
66 #include <Geom2d_BezierCurve.hxx>
67 #include <Geom2d_TrimmedCurve.hxx>
68 #include <GeomAPI_ProjectPointOnSurf.hxx>
69 #include <GeomLib.hxx>
70 #include <Geom_BSplineCurve.hxx>
71 #include <Geom_BSplineSurface.hxx>
72 #include <Geom_BezierCurve.hxx>
73 #include <Geom_BezierSurface.hxx>
74 #include <Geom_RectangularTrimmedSurface.hxx>
75 #include <Geom_TrimmedCurve.hxx>
76 #include <IntAna_IntConicQuad.hxx>
77 #include <IntAna_IntLinTorus.hxx>
78 #include <IntAna_Quadric.hxx>
79 #include <IntCurveSurface_TransitionOnCurve.hxx>
80 #include <IntCurvesFace_Intersector.hxx>
81 #include <Poly_Triangulation.hxx>
82 #include <Precision.hxx>
83 #include <TopExp.hxx>
84 #include <TopExp_Explorer.hxx>
85 #include <TopLoc_Location.hxx>
86 #include <TopTools_DataMapOfShapeInteger.hxx>
87 #include <TopTools_IndexedMapOfShape.hxx>
88 #include <TopTools_MapOfShape.hxx>
89 #include <TopoDS.hxx>
90 #include <TopoDS_Compound.hxx>
91 #include <TopoDS_Face.hxx>
92 #include <TopoDS_TShape.hxx>
93 #include <gp_Cone.hxx>
94 #include <gp_Cylinder.hxx>
95 #include <gp_Lin.hxx>
96 #include <gp_Pln.hxx>
97 #include <gp_Pnt2d.hxx>
98 #include <gp_Sphere.hxx>
99 #include <gp_Torus.hxx>
100
101 //STD
102 #include <limits>
103 #include <mutex>
104 #include <thread>
105
106 #include <boost/container/flat_map.hpp>
107
108 #ifdef _DEBUG_
109 //  #define _MY_DEBUG_
110 //  #undef WITH_TBB
111 #endif
112
113 #ifdef WITH_TBB
114
115 #ifdef WIN32
116 // See https://docs.microsoft.com/en-gb/cpp/porting/modifying-winver-and-win32-winnt?view=vs-2019
117 // Windows 10 = 0x0A00  
118 #define WINVER 0x0A00
119 #define _WIN32_WINNT 0x0A00
120 #endif
121 #include <algorithm>
122 #include <tbb/parallel_for.h>
123 #endif
124
125 using namespace std;
126 using namespace SMESH;
127 std::mutex _eMutex;
128 std::mutex _bMutex;
129
130 //=============================================================================
131 /*!
132  * Constructor
133  */
134 //=============================================================================
135
136 StdMeshers_Cartesian_3D::StdMeshers_Cartesian_3D(int hypId, SMESH_Gen * gen)
137   :SMESH_3D_Algo(hypId, gen)
138 {
139   _name = "Cartesian_3D";
140   _shapeType = (1 << TopAbs_SOLID);       // 1 bit /shape type
141   _compatibleHypothesis.push_back( "CartesianParameters3D" );
142   _compatibleHypothesis.push_back( StdMeshers_ViscousLayers::GetHypType() );
143
144   _onlyUnaryInput = false;          // to mesh all SOLIDs at once
145   _requireDiscreteBoundary = false; // 2D mesh not needed
146   _supportSubmeshes = false;        // do not use any existing mesh
147 }
148
149 //=============================================================================
150 /*!
151  * Check presence of a hypothesis
152  */
153 //=============================================================================
154
155 bool StdMeshers_Cartesian_3D::CheckHypothesis (SMESH_Mesh&          aMesh,
156                                                const TopoDS_Shape&  aShape,
157                                                Hypothesis_Status&   aStatus)
158 {
159   aStatus = SMESH_Hypothesis::HYP_MISSING;
160
161   const list<const SMESHDS_Hypothesis*>& hyps = GetUsedHypothesis(aMesh, aShape, /*skipAux=*/false);
162   list <const SMESHDS_Hypothesis* >::const_iterator h = hyps.begin();
163   if ( h == hyps.end())
164   {
165     return false;
166   }
167
168   _hyp = nullptr;
169   _hypViscousLayers = nullptr;
170   _isComputeOffset = false;
171
172   for ( ; h != hyps.end(); ++h )
173   {
174     if ( !_hyp && ( _hyp = dynamic_cast<const StdMeshers_CartesianParameters3D*>( *h )))
175     {
176       aStatus = _hyp->IsDefined() ? HYP_OK : HYP_BAD_PARAMETER;
177     }
178     else
179     {
180       _hypViscousLayers = dynamic_cast<const StdMeshers_ViscousLayers*>( *h );
181     }
182   }
183
184   return aStatus == HYP_OK;
185 }
186
187 namespace
188 {
189   /*!
190    * \brief Temporary mesh to hold 
191    */
192   struct TmpMesh: public SMESH_Mesh
193   {
194     TmpMesh() {
195       _isShapeToMesh = (_id = 0);
196       _meshDS  = new SMESHDS_Mesh( _id, true );
197     }
198   };
199
200   typedef int                     TGeomID; // IDs of sub-shapes
201   typedef TopTools_ShapeMapHasher TShapeHasher; // non-oriented shape hasher
202   typedef std::array< int, 3 >    TIJK;
203
204   const TGeomID theUndefID = 1e+9;
205
206   //=============================================================================
207   // Definitions of internal utils
208   // --------------------------------------------------------------------------
209   enum Transition {
210     Trans_TANGENT = IntCurveSurface_Tangent,
211     Trans_IN      = IntCurveSurface_In,
212     Trans_OUT     = IntCurveSurface_Out,
213     Trans_APEX,
214     Trans_INTERNAL // for INTERNAL FACE
215   };
216   // --------------------------------------------------------------------------
217   /*!
218    * \brief Sub-entities of a FACE neighboring its concave VERTEX.
219    *        Help to avoid linking nodes on EDGEs that seem connected
220    *        by the concave FACE but the link actually lies outside the FACE
221    */
222   struct ConcaveFace
223   {
224     TGeomID _concaveFace;
225     TGeomID _edge1, _edge2;
226     TGeomID _v1,    _v2;
227     ConcaveFace( int f=0, int e1=0, int e2=0, int v1=0, int v2=0 )
228       : _concaveFace(f), _edge1(e1), _edge2(e2), _v1(v1), _v2(v2) {}
229     bool HasEdge( TGeomID edge ) const { return edge == _edge1 || edge == _edge2; }
230     bool HasVertex( TGeomID v  ) const { return v == _v1 || v == _v2; }
231     void SetEdge( TGeomID edge ) { ( _edge1 ? _edge2 : _edge1 ) = edge; }
232     void SetVertex( TGeomID v  ) { ( _v1 ? _v2 : _v1 ) = v; }
233   };
234   typedef NCollection_DataMap< TGeomID, ConcaveFace > TConcaveVertex2Face;
235   // --------------------------------------------------------------------------
236   /*!
237    * \brief Container of IDs of SOLID sub-shapes
238    */
239   class Solid // sole SOLID contains all sub-shapes
240   {
241     TGeomID             _id; // SOLID id
242     bool                _hasInternalFaces;
243     TConcaveVertex2Face _concaveVertex; // concave VERTEX -> ConcaveFace
244   public:
245     virtual ~Solid() {}
246     virtual bool Contains( TGeomID /*subID*/ ) const { return true; }
247     virtual bool ContainsAny( const vector< TGeomID>& /*subIDs*/ ) const { return true; }
248     virtual TopAbs_Orientation Orientation( const TopoDS_Shape& s ) const { return s.Orientation(); }
249     virtual bool IsOutsideOriented( TGeomID /*faceID*/ ) const { return true; }
250     void SetID( TGeomID id ) { _id = id; }
251     TGeomID ID() const { return _id; }
252     void SetHasInternalFaces( bool has ) { _hasInternalFaces = has; }
253     bool HasInternalFaces() const { return _hasInternalFaces; }
254     void SetConcave( TGeomID V, TGeomID F, TGeomID E1, TGeomID E2, TGeomID V1, TGeomID V2  )
255     { _concaveVertex.Bind( V, ConcaveFace{ F, E1, E2, V1, V2 }); }
256     bool HasConcaveVertex() const { return !_concaveVertex.IsEmpty(); }
257     const ConcaveFace* GetConcave( TGeomID V ) const { return _concaveVertex.Seek( V ); }
258   };
259   // --------------------------------------------------------------------------
260   class OneOfSolids : public Solid
261   {
262     TColStd_MapOfInteger _subIDs;
263     TopTools_MapOfShape  _faces; // keep FACE orientation
264     TColStd_MapOfInteger _outFaceIDs; // FACEs of shape_to_mesh oriented outside the SOLID
265   public:
266     void Init( const TopoDS_Shape& solid,
267                TopAbs_ShapeEnum    subType,
268                const SMESHDS_Mesh* mesh );
269     virtual bool Contains( TGeomID i ) const { return i == ID() || _subIDs.Contains( i ); }
270     virtual bool ContainsAny( const vector< TGeomID>& subIDs ) const
271     {
272       for ( size_t i = 0; i < subIDs.size(); ++i ) if ( Contains( subIDs[ i ])) return true;
273       return false;
274     }
275     virtual TopAbs_Orientation Orientation( const TopoDS_Shape& face ) const
276     {
277       const TopoDS_Shape& sInMap = const_cast< OneOfSolids* >(this)->_faces.Added( face );
278       return sInMap.Orientation();
279     }
280     virtual bool IsOutsideOriented( TGeomID faceID ) const
281     {
282       return faceID == 0 || _outFaceIDs.Contains( faceID );
283     }
284   };
285   // --------------------------------------------------------------------------
286   /*!
287    * \brief Hold a vector of TGeomID and clear it at destruction
288    */
289   class GeomIDVecHelder
290   {
291     typedef std::vector< TGeomID > TVector;
292     const TVector& myVec;
293     bool           myOwn;
294
295   public:
296     GeomIDVecHelder( const TVector& idVec, bool isOwner ): myVec( idVec ), myOwn( isOwner ) {}
297     GeomIDVecHelder( const GeomIDVecHelder& holder ): myVec( holder.myVec ), myOwn( holder.myOwn )
298     {
299       const_cast< bool& >( holder.myOwn ) = false;
300     }
301     ~GeomIDVecHelder() { if ( myOwn ) const_cast<TVector&>( myVec ).clear(); }
302     size_t size() const { return myVec.size(); }
303     TGeomID operator[]( size_t i ) const { return i < size() ? myVec[i] : theUndefID; }
304     bool operator==( const GeomIDVecHelder& other ) const { return myVec == other.myVec; }
305     bool contain( const TGeomID& id ) const {
306       return std::find( myVec.begin(), myVec.end(), id ) != myVec.end();
307     }
308     TGeomID otherThan( const TGeomID& id ) const {
309       for ( const TGeomID& id2 : myVec )
310         if ( id != id2 )
311           return id2;
312       return theUndefID;
313     }
314     TGeomID oneCommon( const GeomIDVecHelder& other ) const {
315       TGeomID common = theUndefID;
316       for ( const TGeomID& id : myVec )
317         if ( other.contain( id ))
318         {
319           if ( common != theUndefID )
320             return theUndefID;
321           common = id;
322         }
323       return common;
324     }
325   };
326   // --------------------------------------------------------------------------
327   /*!
328    * \brief Geom data
329    */
330   struct Geometry
331   {
332     TopoDS_Shape                _mainShape;
333     vector< vector< TGeomID > > _solidIDsByShapeID;// V/E/F ID -> SOLID IDs
334     Solid                       _soleSolid;
335     map< TGeomID, OneOfSolids > _solidByID;
336     TColStd_MapOfInteger        _boundaryFaces; // FACEs on boundary of mesh->ShapeToMesh()
337     TColStd_MapOfInteger        _strangeEdges; // EDGEs shared by strange FACEs
338     TGeomID                     _extIntFaceID; // pseudo FACE - extension of INTERNAL FACE
339
340     TopTools_DataMapOfShapeInteger _shape2NbNodes; // nb of pre-existing nodes on shapes
341
342     Controls::ElementsOnShape _edgeClassifier;
343     Controls::ElementsOnShape _vertexClassifier;
344
345     bool IsOneSolid() const { return _solidByID.size() < 2; }
346     GeomIDVecHelder GetSolidIDsByShapeID( const vector< TGeomID >& shapeIDs ) const;
347   };
348   // --------------------------------------------------------------------------
349   /*!
350    * \brief Common data of any intersection between a Grid and a shape
351    */
352   struct B_IntersectPoint
353   {
354     mutable const SMDS_MeshNode* _node;
355     mutable vector< TGeomID >    _faceIDs;
356
357     B_IntersectPoint(): _node(NULL) {}
358     bool Add( const vector< TGeomID >& fIDs, const SMDS_MeshNode* n=NULL ) const;
359     TGeomID HasCommonFace( const B_IntersectPoint * other, TGeomID avoidFace=-1 ) const;
360     size_t GetCommonFaces( const B_IntersectPoint * other, TGeomID * commonFaces ) const;
361     bool IsOnFace( TGeomID faceID ) const;
362     virtual ~B_IntersectPoint() {}
363   };
364   // --------------------------------------------------------------------------
365   /*!
366    * \brief Data of intersection between a GridLine and a TopoDS_Face
367    */
368   struct F_IntersectPoint : public B_IntersectPoint
369   {
370     double             _paramOnLine;
371     double             _u, _v;
372     mutable Transition _transition;
373     mutable size_t     _indexOnLine;
374
375     bool operator< ( const F_IntersectPoint& o ) const { return _paramOnLine < o._paramOnLine; }
376   };
377   // --------------------------------------------------------------------------
378   /*!
379    * \brief Data of intersection between GridPlanes and a TopoDS_EDGE
380    */
381   struct E_IntersectPoint : public B_IntersectPoint
382   {
383     gp_Pnt  _point;
384     double  _uvw[3];
385     TGeomID _shapeID; // ID of EDGE or VERTEX
386   };
387   // --------------------------------------------------------------------------
388   /*!
389    * \brief A line of the grid and its intersections with 2D geometry
390    */
391   struct GridLine
392   {
393     gp_Lin _line;
394     double _length; // line length
395     multiset< F_IntersectPoint > _intPoints;
396
397     void RemoveExcessIntPoints( const double tol );
398     TGeomID GetSolidIDBefore( multiset< F_IntersectPoint >::iterator ip,
399                               const TGeomID                          prevID,
400                               const Geometry&                        geom);
401   };
402   // --------------------------------------------------------------------------
403   /*!
404    * \brief Planes of the grid used to find intersections of an EDGE with a hexahedron
405    */
406   struct GridPlanes
407   {
408     gp_XYZ           _zNorm;
409     vector< gp_XYZ > _origins; // origin points of all planes in one direction
410     vector< double > _zProjs;  // projections of origins to _zNorm
411   };
412   // --------------------------------------------------------------------------
413   /*!
414    * \brief Iterator on the parallel grid lines of one direction
415    */
416   struct LineIndexer
417   {
418     size_t _size  [3];
419     size_t _curInd[3];
420     size_t _iVar1, _iVar2, _iConst;
421     string _name1, _name2, _nameConst;
422     LineIndexer() {}
423     LineIndexer( size_t sz1, size_t sz2, size_t sz3,
424                  size_t iv1, size_t iv2, size_t iConst,
425                  const string& nv1, const string& nv2, const string& nConst )
426     {
427       _size[0] = sz1; _size[1] = sz2; _size[2] = sz3;
428       _curInd[0] = _curInd[1] = _curInd[2] = 0;
429       _iVar1 = iv1; _iVar2 = iv2; _iConst = iConst;
430       _name1 = nv1; _name2 = nv2; _nameConst = nConst;
431     }
432
433     size_t I() const { return _curInd[0]; }
434     size_t J() const { return _curInd[1]; }
435     size_t K() const { return _curInd[2]; }
436     void SetIJK( size_t i, size_t j, size_t k )
437     {
438       _curInd[0] = i; _curInd[1] = j; _curInd[2] = k;
439     }
440     void SetLineIndex(size_t i)
441     {
442       _curInd[_iVar2] = i / _size[_iVar1];
443       _curInd[_iVar1] = i % _size[_iVar1];
444     }
445     void operator++()
446     {
447       if ( ++_curInd[_iVar1] == _size[_iVar1] )
448         _curInd[_iVar1] = 0, ++_curInd[_iVar2];
449     }
450     bool More() const { return _curInd[_iVar2] < _size[_iVar2]; }
451     size_t LineIndex   () const { return _curInd[_iVar1] + _curInd[_iVar2]* _size[_iVar1]; }
452     size_t LineIndex10 () const { return (_curInd[_iVar1] + 1 ) + _curInd[_iVar2]* _size[_iVar1]; }
453     size_t LineIndex01 () const { return _curInd[_iVar1] + (_curInd[_iVar2] + 1 )* _size[_iVar1]; }
454     size_t LineIndex11 () const { return (_curInd[_iVar1] + 1 ) + (_curInd[_iVar2] + 1 )* _size[_iVar1]; }
455     void SetIndexOnLine (size_t i)  { _curInd[ _iConst ] = i; }
456     bool IsValidIndexOnLine (size_t i) const { return  i < _size[ _iConst ]; }
457     size_t NbLines() const { return _size[_iVar1] * _size[_iVar2]; }
458   };
459   struct FaceGridIntersector;
460   // --------------------------------------------------------------------------
461   /*!
462    * \brief Container of GridLine's
463    */
464   struct Grid
465   {
466     vector< double >   _coords[3]; // coordinates of grid nodes
467     gp_XYZ             _axes  [3]; // axis directions
468     vector< GridLine > _lines [3]; //    in 3 directions
469     double             _tol, _minCellSize;
470     gp_XYZ             _origin;
471     gp_Mat             _invB; // inverted basis of _axes
472
473     // index shift within _nodes of nodes of a cell from the 1st node
474     int                _nodeShift[8];
475
476     vector< const SMDS_MeshNode* >    _nodes;          // mesh nodes at grid nodes
477     vector< const SMDS_MeshNode* >    _allBorderNodes; // mesh nodes between the bounding box and the geometry boundary
478
479     vector< const F_IntersectPoint* > _gridIntP; // grid node intersection with geometry
480     ObjectPool< E_IntersectPoint >    _edgeIntPool; // intersections with EDGEs
481     ObjectPool< F_IntersectPoint >    _extIntPool; // intersections with extended INTERNAL FACEs
482     //list< E_IntersectPoint >          _edgeIntP; // intersections with EDGEs
483
484     Geometry                          _geometry;
485     bool                              _toAddEdges;
486     bool                              _toCreateFaces;
487     bool                              _toConsiderInternalFaces;
488     bool                              _toUseThresholdForInternalFaces;
489     double                            _sizeThreshold;
490     bool                              _toUseQuanta;
491     double                            _quanta;
492
493     SMESH_MesherHelper*               _helper;
494
495     size_t CellIndex( size_t i, size_t j, size_t k ) const
496     {
497       return i + j*(_coords[0].size()-1) + k*(_coords[0].size()-1)*(_coords[1].size()-1);
498     }
499     size_t NodeIndex( size_t i, size_t j, size_t k ) const
500     {
501       return i + j*_coords[0].size() + k*_coords[0].size()*_coords[1].size();
502     }
503     size_t NodeIndex( const TIJK& ijk ) const
504     {
505       return NodeIndex( ijk[0], ijk[1], ijk[2] );
506     }
507     size_t NodeIndexDX() const { return 1; }
508     size_t NodeIndexDY() const { return _coords[0].size(); }
509     size_t NodeIndexDZ() const { return _coords[0].size() * _coords[1].size(); }
510
511     LineIndexer GetLineIndexer(size_t iDir) const;
512     size_t GetLineDir( const GridLine* line, size_t & index ) const;
513
514     E_IntersectPoint* Add( const E_IntersectPoint& ip )
515     {
516       E_IntersectPoint* eip = _edgeIntPool.getNew();
517       *eip = ip;
518       return eip;
519     }
520     void Remove( E_IntersectPoint* eip ) { _edgeIntPool.destroy( eip ); }
521
522     TGeomID ShapeID( const TopoDS_Shape& s ) const;
523     const TopoDS_Shape& Shape( TGeomID id ) const;
524     TopAbs_ShapeEnum ShapeType( TGeomID id ) const { return Shape(id).ShapeType(); }
525     void InitGeometry( const TopoDS_Shape& theShape );
526     void InitClassifier( const TopoDS_Shape&        mainShape,
527                          TopAbs_ShapeEnum           shapeType,
528                          Controls::ElementsOnShape& classifier );
529     void GetEdgesToImplement( map< TGeomID, vector< TGeomID > > & edge2faceMap,
530                               const TopoDS_Shape&                 shape,
531                               const vector< TopoDS_Shape >&       faces );
532     void SetSolidFather( const TopoDS_Shape& s, const TopoDS_Shape& theShapeToMesh );
533     bool IsShared( TGeomID faceID ) const;
534     bool IsAnyShared( const std::vector< TGeomID >& faceIDs ) const;
535     bool IsInternal( TGeomID faceID ) const {
536       return ( faceID == PseudoIntExtFaceID() ||
537                Shape( faceID ).Orientation() == TopAbs_INTERNAL ); }
538     bool IsSolid( TGeomID shapeID ) const {
539       if ( _geometry.IsOneSolid() ) return _geometry._soleSolid.ID() == shapeID;
540       else                          return _geometry._solidByID.count( shapeID ); }
541     bool IsStrangeEdge( TGeomID id ) const { return _geometry._strangeEdges.Contains( id ); }
542     TGeomID PseudoIntExtFaceID() const { return _geometry._extIntFaceID; }
543     Solid* GetSolid( TGeomID solidID = 0 );
544     Solid* GetOneOfSolids( TGeomID solidID );
545     const vector< TGeomID > & GetSolidIDs( TGeomID subShapeID ) const;
546     bool IsCorrectTransition( TGeomID faceID, const Solid* solid );
547     bool IsBoundaryFace( TGeomID face ) const { return _geometry._boundaryFaces.Contains( face ); }
548     void SetOnShape( const SMDS_MeshNode* n, const F_IntersectPoint& ip,
549                      TopoDS_Vertex* vertex = nullptr, bool unset = false );
550     void UpdateFacesOfVertex( const B_IntersectPoint& ip, const TopoDS_Vertex& vertex );
551     bool IsToCheckNodePos() const { return !_toAddEdges && _toCreateFaces; }
552     bool IsToRemoveExcessEntities() const { return !_toAddEdges; }
553
554     void SetCoordinates(const vector<double>& xCoords,
555                         const vector<double>& yCoords,
556                         const vector<double>& zCoords,
557                         const double*         axesDirs,
558                         const Bnd_Box&        bndBox );
559     void ComputeUVW(const gp_XYZ& p, double uvw[3]);
560     void ComputeNodes(SMESH_MesherHelper& helper);
561   };
562   // --------------------------------------------------------------------------
563   /*!
564    * \brief Return cells sharing a link
565    */
566   struct CellsAroundLink
567   {
568     int    _iDir;
569     int    _dInd[4][3];
570     size_t _nbCells[3];
571     int    _i,_j,_k;
572     Grid*  _grid;
573
574     CellsAroundLink( Grid* grid, int iDir ):
575       _iDir( iDir ),
576       _dInd{ {0,0,0}, {0,0,0}, {0,0,0}, {0,0,0} },
577       _nbCells{ grid->_coords[0].size() - 1,
578           grid->_coords[1].size() - 1,
579           grid->_coords[2].size() - 1 },
580       _grid( grid )
581     {
582       const int iDirOther[3][2] = {{ 1,2 },{ 0,2 },{ 0,1 }};
583       _dInd[1][ iDirOther[iDir][0] ] = -1;
584       _dInd[2][ iDirOther[iDir][1] ] = -1;
585       _dInd[3][ iDirOther[iDir][0] ] = -1; _dInd[3][ iDirOther[iDir][1] ] = -1;
586     }
587     void Init( int i, int j, int k, int link12 = 0 )
588     {
589       int iL = link12 % 4;
590       _i = i - _dInd[iL][0];
591       _j = j - _dInd[iL][1];
592       _k = k - _dInd[iL][2];
593     }
594     bool GetCell( int iL, int& i, int& j, int& k, int& cellIndex, int& linkIndex )
595     {
596       i =  _i + _dInd[iL][0];
597       j =  _j + _dInd[iL][1];
598       k =  _k + _dInd[iL][2];
599       if ( i < 0 || i >= (int)_nbCells[0] ||
600            j < 0 || j >= (int)_nbCells[1] ||
601            k < 0 || k >= (int)_nbCells[2] )
602         return false;
603       cellIndex = _grid->CellIndex( i,j,k );
604       linkIndex = iL + _iDir * 4;
605       return true;
606     }
607   };
608   // --------------------------------------------------------------------------
609   /*!
610    * \brief Intersector of TopoDS_Face with all GridLine's
611    */
612   struct FaceGridIntersector
613   {
614     TopoDS_Face _face;
615     TGeomID     _faceID;
616     Grid*       _grid;
617     Bnd_Box     _bndBox;
618     IntCurvesFace_Intersector* _surfaceInt;
619     vector< std::pair< GridLine*, F_IntersectPoint > > _intersections;
620
621     FaceGridIntersector(): _grid(0), _surfaceInt(0) {}
622     void Intersect();
623
624     void StoreIntersections()
625     {
626       for ( size_t i = 0; i < _intersections.size(); ++i )
627       {
628         multiset< F_IntersectPoint >::iterator ip =
629           _intersections[i].first->_intPoints.insert( _intersections[i].second );
630         ip->_faceIDs.reserve( 1 );
631         ip->_faceIDs.push_back( _faceID );
632       }
633     }
634     const Bnd_Box& GetFaceBndBox()
635     {
636       GetCurveFaceIntersector();
637       return _bndBox;
638     }
639     IntCurvesFace_Intersector* GetCurveFaceIntersector()
640     {
641       if ( !_surfaceInt )
642       {
643         _surfaceInt = new IntCurvesFace_Intersector( _face, Precision::PConfusion() );
644         _bndBox     = _surfaceInt->Bounding();
645         if ( _bndBox.IsVoid() )
646           BRepBndLib::Add (_face, _bndBox);
647       }
648       return _surfaceInt;
649     }
650     bool IsThreadSafe(set< const Standard_Transient* >& noSafeTShapes) const;
651   };
652   // --------------------------------------------------------------------------
653   /*!
654    * \brief Intersector of a surface with a GridLine
655    */
656   struct FaceLineIntersector
657   {
658     double      _tol;
659     double      _u, _v, _w; // params on the face and the line
660     Transition  _transition; // transition at intersection (see IntCurveSurface.cdl)
661     Transition  _transIn, _transOut; // IN and OUT transitions depending of face orientation
662
663     gp_Pln      _plane;
664     gp_Cylinder _cylinder;
665     gp_Cone     _cone;
666     gp_Sphere   _sphere;
667     gp_Torus    _torus;
668     IntCurvesFace_Intersector* _surfaceInt;
669
670     vector< F_IntersectPoint > _intPoints;
671
672     void IntersectWithPlane   (const GridLine& gridLine);
673     void IntersectWithCylinder(const GridLine& gridLine);
674     void IntersectWithCone    (const GridLine& gridLine);
675     void IntersectWithSphere  (const GridLine& gridLine);
676     void IntersectWithTorus   (const GridLine& gridLine);
677     void IntersectWithSurface (const GridLine& gridLine);
678
679     bool UVIsOnFace() const;
680     void addIntPoint(const bool toClassify=true);
681     bool isParamOnLineOK( const double linLength )
682     {
683       return -_tol < _w && _w < linLength + _tol;
684     }
685     FaceLineIntersector():_surfaceInt(0) {}
686     ~FaceLineIntersector() { if (_surfaceInt ) delete _surfaceInt; _surfaceInt = 0; }
687   };
688   // --------------------------------------------------------------------------
689   /*!
690    * \brief Class representing topology of the hexahedron and creating a mesh
691    *        volume basing on analysis of hexahedron intersection with geometry
692    */
693   class Hexahedron
694   {
695     // --------------------------------------------------------------------------------
696     struct _Face;
697     struct _Link;
698     enum IsInternalFlag { IS_NOT_INTERNAL, IS_INTERNAL, IS_CUT_BY_INTERNAL_FACE };
699     // --------------------------------------------------------------------------------
700     struct _Node //!< node either at a hexahedron corner or at intersection
701     {
702       const SMDS_MeshNode*    _node;        // mesh node at hexahedron corner
703       const SMDS_MeshNode*    _boundaryCornerNode; // missing mesh node due to hex truncation on the boundary
704       const B_IntersectPoint* _intPoint;
705       const _Face*            _usedInFace;
706       char                    _isInternalFlags;
707
708       _Node(const SMDS_MeshNode* n=0, const B_IntersectPoint* ip=0)
709         :_node(n), _intPoint(ip), _usedInFace(0), _isInternalFlags(0) {} 
710       const SMDS_MeshNode*    Node() const
711       { return ( _intPoint && _intPoint->_node ) ? _intPoint->_node : _node; }
712       const SMDS_MeshNode*    BoundaryNode() const
713       { return _node ? _node : _boundaryCornerNode; }
714       const E_IntersectPoint* EdgeIntPnt() const
715       { return static_cast< const E_IntersectPoint* >( _intPoint ); }
716       const F_IntersectPoint* FaceIntPnt() const
717       { return static_cast< const F_IntersectPoint* >( _intPoint ); }
718       const vector< TGeomID >& faces() const { return _intPoint->_faceIDs; }
719       TGeomID face(size_t i) const { return _intPoint->_faceIDs[ i ]; }
720       void SetInternal( IsInternalFlag intFlag ) { _isInternalFlags |= intFlag; }
721       bool IsCutByInternal() const { return _isInternalFlags & IS_CUT_BY_INTERNAL_FACE; }
722       bool IsUsedInFace( const _Face* polygon = 0 )
723       {
724         return polygon ? ( _usedInFace == polygon ) : bool( _usedInFace );
725       }
726       TGeomID IsLinked( const B_IntersectPoint* other,
727                         TGeomID                 avoidFace=-1 ) const // returns id of a common face
728       {
729         return _intPoint ? _intPoint->HasCommonFace( other, avoidFace ) : 0;
730       }
731       bool IsOnFace( TGeomID faceID ) const // returns true if faceID is found
732       {
733         return _intPoint ? _intPoint->IsOnFace( faceID ) : false;
734       }
735       size_t GetCommonFaces( const B_IntersectPoint * other, TGeomID* common ) const
736       {
737         return _intPoint && other ? _intPoint->GetCommonFaces( other, common ) : 0;
738       }
739       gp_Pnt Point() const
740       {
741         if ( const SMDS_MeshNode* n = Node() )
742           return SMESH_NodeXYZ( n );
743         if ( const E_IntersectPoint* eip =
744              dynamic_cast< const E_IntersectPoint* >( _intPoint ))
745           return eip->_point;
746         return gp_Pnt( 1e100, 0, 0 );
747       }
748       TGeomID ShapeID() const
749       {
750         if ( const E_IntersectPoint* eip = dynamic_cast< const E_IntersectPoint* >( _intPoint ))
751           return eip->_shapeID;
752         return 0;
753       }
754       void Add( const E_IntersectPoint* ip )
755       {
756         const std::lock_guard<std::mutex> lock(_eMutex);
757         // Possible cases before Add(ip):
758         ///  1) _node != 0 --> _Node at hex corner ( _intPoint == 0 || _intPoint._node == 0 )
759         ///  2) _node == 0 && _intPoint._node != 0  -->  link intersected by FACE
760         ///  3) _node == 0 && _intPoint._node == 0  -->  _Node at EDGE intersection
761         //
762         // If ip is added in cases 1) and 2) _node position must be changed to ip._shapeID
763         //   at creation of elements
764         // To recognize this case, set _intPoint._node = Node()
765         const SMDS_MeshNode* node = Node();
766         if ( !_intPoint ) {
767           _intPoint = ip;
768         }
769         else {
770           ip->Add( _intPoint->_faceIDs );
771           _intPoint = ip;
772         }
773         if ( node )
774           _node = _intPoint->_node = node;
775       }
776     };
777     // --------------------------------------------------------------------------------
778     struct _Link // link connecting two _Node's
779     {
780       _Node* _nodes[2];
781       _Face* _faces[2]; // polygons sharing a link
782       vector< const F_IntersectPoint* > _fIntPoints; // GridLine intersections with FACEs
783       vector< _Node* >                  _fIntNodes;   // _Node's at _fIntPoints
784       vector< _Link >                   _splits;
785       _Link(): _faces{ 0, 0 } {}
786     };
787     // --------------------------------------------------------------------------------
788     struct _OrientedLink
789     {
790       _Link* _link;
791       bool   _reverse;
792       _OrientedLink( _Link* link=0, bool reverse=false ): _link(link), _reverse(reverse) {}
793       void Reverse() { _reverse = !_reverse; }
794       size_t NbResultLinks() const { return _link->_splits.size(); }
795       _OrientedLink ResultLink(int i) const
796       {
797         return _OrientedLink(&_link->_splits[_reverse ? NbResultLinks()-i-1 : i],_reverse);
798       }
799       _Node* FirstNode() const { return _link->_nodes[ _reverse ]; }
800       _Node* LastNode()  const { return _link->_nodes[ !_reverse ]; }
801       operator bool() const { return _link; }
802       vector< TGeomID > GetNotUsedFace(const set<TGeomID>& usedIDs ) const // returns supporting FACEs
803       {
804         vector< TGeomID > faces;
805         const B_IntersectPoint *ip0, *ip1;
806         if (( ip0 = _link->_nodes[0]->_intPoint ) &&
807             ( ip1 = _link->_nodes[1]->_intPoint ))
808         {
809           for ( size_t i = 0; i < ip0->_faceIDs.size(); ++i )
810             if ( ip1->IsOnFace ( ip0->_faceIDs[i] ) &&
811                  !usedIDs.count( ip0->_faceIDs[i] ) )
812               faces.push_back( ip0->_faceIDs[i] );
813         }
814         return faces;
815       }
816       bool HasEdgeNodes() const
817       {
818         return ( dynamic_cast< const E_IntersectPoint* >( _link->_nodes[0]->_intPoint ) ||
819                  dynamic_cast< const E_IntersectPoint* >( _link->_nodes[1]->_intPoint ));
820       }
821       int NbFaces() const
822       {
823         return !_link->_faces[0] ? 0 : 1 + bool( _link->_faces[1] );
824       }
825       void AddFace( _Face* f )
826       {
827         if ( _link->_faces[0] )
828         {
829           _link->_faces[1] = f;
830         }
831         else
832         {
833           _link->_faces[0] = f;
834           _link->_faces[1] = 0;
835         }
836       }
837       void RemoveFace( _Face* f )
838       {
839         if ( !_link->_faces[0] ) return;
840
841         if ( _link->_faces[1] == f )
842         {
843           _link->_faces[1] = 0;
844         }
845         else if ( _link->_faces[0] == f )
846         {
847           _link->_faces[0] = 0;
848           if ( _link->_faces[1] )
849           {
850             _link->_faces[0] = _link->_faces[1];
851             _link->_faces[1] = 0;
852           }
853         }
854       }
855     };
856     // --------------------------------------------------------------------------------
857     struct _SplitIterator //! set to _hexLinks splits on one side of INTERNAL FACEs
858     {
859       struct _Split // data of a link split
860       {
861         int    _linkID;          // hex link ID
862         _Node* _nodes[2];
863         int    _iCheckIteration; // iteration where split is tried as Hexahedron split
864         _Link* _checkedSplit;    // split set to hex links
865         bool   _isUsed;          // used in a volume
866
867         _Split( _Link & split, int iLink ):
868           _linkID( iLink ), _nodes{ split._nodes[0], split._nodes[1] },
869           _iCheckIteration( 0 ), _isUsed( false )
870         {}
871         bool IsCheckedOrUsed( bool used ) const { return used ? _isUsed : _iCheckIteration > 0; }
872       };
873       _Link*                _hexLinks;
874       std::vector< _Split > _splits;
875       int                   _iterationNb;
876       size_t                _nbChecked;
877       size_t                _nbUsed;
878       std::vector< _Node* > _freeNodes; // nodes reached while composing a split set
879
880       _SplitIterator( _Link* hexLinks ):
881         _hexLinks( hexLinks ), _iterationNb(0), _nbChecked(0), _nbUsed(0)
882       {
883         _freeNodes.reserve( 12 );
884         _splits.reserve( 24 );
885         for ( int iL = 0; iL < 12; ++iL )
886           for ( size_t iS = 0; iS < _hexLinks[ iL ]._splits.size(); ++iS )
887             _splits.emplace_back( _hexLinks[ iL ]._splits[ iS ], iL );
888         Next();
889       }
890       bool More() const { return _nbUsed < _splits.size(); }
891       bool Next();
892     };
893     // --------------------------------------------------------------------------------
894     struct _Face
895     {
896       SMESH_Block::TShapeID   _name;
897       vector< _OrientedLink > _links;       // links on GridLine's
898       vector< _Link >         _polyLinks;   // links added to close a polygonal face
899       vector< _Node* >        _eIntNodes;   // nodes at intersection with EDGEs
900
901       _Face():_name( SMESH_Block::ID_NONE )
902       {}
903       bool IsPolyLink( const _OrientedLink& ol )
904       {
905         return _polyLinks.empty() ? false :
906           ( &_polyLinks[0] <= ol._link &&  ol._link <= &_polyLinks.back() );
907       }
908       void AddPolyLink(_Node* n0, _Node* n1, _Face* faceToFindEqual=0)
909       {
910         if ( faceToFindEqual && faceToFindEqual != this ) {
911           for ( size_t iL = 0; iL < faceToFindEqual->_polyLinks.size(); ++iL )
912             if ( faceToFindEqual->_polyLinks[iL]._nodes[0] == n1 &&
913                  faceToFindEqual->_polyLinks[iL]._nodes[1] == n0 )
914             {
915               _links.push_back
916                 ( _OrientedLink( & faceToFindEqual->_polyLinks[iL], /*reverse=*/true ));
917               return;
918             }
919         }
920         _Link l;
921         l._nodes[0] = n0;
922         l._nodes[1] = n1;
923         _polyLinks.push_back( l );
924         _links.push_back( _OrientedLink( &_polyLinks.back() ));
925       }
926     };
927     // --------------------------------------------------------------------------------
928     struct _volumeDef // holder of nodes of a volume mesh element
929     {
930       typedef void* _ptr;
931
932       struct _nodeDef
933       {
934         const SMDS_MeshNode*    _node; // mesh node at hexahedron corner
935         const B_IntersectPoint* _intPoint;
936
937         _nodeDef(): _node(0), _intPoint(0) {}
938         _nodeDef( _Node* n ): _node( n->_node), _intPoint( n->_intPoint ) {}
939         const SMDS_MeshNode*    Node() const
940         { return ( _intPoint && _intPoint->_node ) ? _intPoint->_node : _node; }
941         const E_IntersectPoint* EdgeIntPnt() const
942         { return static_cast< const E_IntersectPoint* >( _intPoint ); }
943         _ptr Ptr() const { return Node() ? (_ptr) Node() : (_ptr) EdgeIntPnt(); }
944         bool operator==(const _nodeDef& other ) const { return Ptr() == other.Ptr(); }
945       };
946
947       vector< _nodeDef >      _nodes;
948       vector< int >           _quantities;
949       _volumeDef*             _next; // to store several _volumeDefs in a chain
950       TGeomID                 _solidID;
951       double                  _size;
952       const SMDS_MeshElement* _volume; // new volume
953       std::vector<const SMDS_MeshElement*> _brotherVolume; // produced due to poly split 
954
955       vector< SMESH_Block::TShapeID > _names; // name of side a polygon originates from
956
957       _volumeDef(): _next(0), _solidID(0), _size(0), _volume(0) {}
958       ~_volumeDef() { delete _next; }
959       _volumeDef( _volumeDef& other ):
960         _next(0), _solidID( other._solidID ), _size( other._size ), _volume( other._volume )
961       { _nodes.swap( other._nodes ); _quantities.swap( other._quantities ); other._volume = 0;
962         _names.swap( other._names ); }
963
964       size_t size() const { return 1 + ( _next ? _next->size() : 0 ); } // nb _volumeDef in a chain
965       _volumeDef* at(int index)
966       { return index == 0 ? this : ( _next ? _next->at(index-1) : _next ); }
967
968       void Set( _Node** nodes, int nb )
969       { _nodes.assign( nodes, nodes + nb ); }
970
971       void SetNext( _volumeDef* vd )
972       { if ( _next ) { _next->SetNext( vd ); } else { _next = vd; }}
973
974       bool IsEmpty() const { return (( _nodes.empty() ) &&
975                                      ( !_next || _next->IsEmpty() )); }
976       bool IsPolyhedron() const { return ( !_quantities.empty() ||
977                                            ( _next && !_next->_quantities.empty() )); }
978
979
980       struct _linkDef: public std::pair<_ptr,_ptr> // to join polygons in removeExcessSideDivision()
981       {
982         _nodeDef _node1;//, _node2;
983         mutable /*const */_linkDef *_prev, *_next;
984         size_t _loopIndex;
985
986         _linkDef():_prev(0), _next(0) {}
987
988         void init( const _nodeDef& n1, const _nodeDef& n2, size_t iLoop )
989         {
990           _node1     = n1; //_node2 = n2;
991           _loopIndex = iLoop;
992           first      = n1.Ptr();
993           second     = n2.Ptr();
994           if ( first > second ) std::swap( first, second );
995         }
996         void setNext( _linkDef* next )
997         {
998           _next = next;
999           next->_prev = this;
1000         }
1001       };
1002     };
1003
1004     // topology of a hexahedron
1005     _Node _hexNodes [8];
1006     _Link _hexLinks [12];
1007     _Face _hexQuads [6];
1008
1009     // faces resulted from hexahedron intersection
1010     vector< _Face > _polygons;
1011
1012     // intresections with EDGEs
1013     vector< const E_IntersectPoint* > _eIntPoints;
1014
1015     // additional nodes created at intersection points
1016     vector< _Node > _intNodes;
1017
1018     // nodes inside the hexahedron (at VERTEXes) refer to _intNodes
1019     vector< _Node* > _vIntNodes;
1020
1021     // computed volume elements
1022     _volumeDef _volumeDefs;
1023
1024     Grid*       _grid;
1025     double      _sideLength[3];
1026     int         _nbCornerNodes, _nbFaceIntNodes, _nbBndNodes;
1027     int         _origNodeInd; // index of _hexNodes[0] node within the _grid
1028     size_t      _i,_j,_k;
1029     bool        _hasTooSmall;
1030     int         _cellID;
1031
1032   public:
1033     Hexahedron(Grid* grid);
1034     int MakeElements(SMESH_MesherHelper&                      helper,
1035                      const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap);
1036     void computeElements( const Solid* solid = 0, int solidIndex = -1 );
1037
1038   private:
1039     Hexahedron(const Hexahedron& other, size_t i, size_t j, size_t k, int cellID );
1040     void init( size_t i, size_t j, size_t k, const Solid* solid=0 );
1041     void init( size_t i );
1042     void setIJK( size_t i );
1043     bool compute( const Solid* solid, const IsInternalFlag intFlag );
1044     size_t getSolids( TGeomID ids[] );
1045     bool isCutByInternalFace( IsInternalFlag & maxFlag );
1046     void addEdges(SMESH_MesherHelper&                      helper,
1047                   vector< Hexahedron* >&                   intersectedHex,
1048                   const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap);
1049     gp_Pnt findIntPoint( double u1, double proj1, double u2, double proj2,
1050                          double proj, BRepAdaptor_Curve& curve,
1051                          const gp_XYZ& axis, const gp_XYZ& origin );
1052     int  getEntity( const E_IntersectPoint* ip, int* facets, int& sub );
1053     bool addIntersection( const E_IntersectPoint* ip,
1054                           vector< Hexahedron* >&  hexes,
1055                           int ijk[], int dIJK[] );
1056     bool isQuadOnFace( const size_t iQuad );
1057     bool findChain( _Node* n1, _Node* n2, _Face& quad, vector<_Node*>& chainNodes );
1058     bool closePolygon( _Face* polygon, vector<_Node*>& chainNodes ) const;
1059     bool findChainOnEdge( const vector< _OrientedLink >& splits,
1060                           const _OrientedLink&           prevSplit,
1061                           const _OrientedLink&           avoidSplit,
1062                           const std::set< TGeomID > &    concaveFaces,
1063                           size_t &                       iS,
1064                           _Face&                         quad,
1065                           vector<_Node*>&                chn);
1066     int  addVolumes(SMESH_MesherHelper& helper );
1067     void addFaces( SMESH_MesherHelper&                       helper,
1068                    const vector< const SMDS_MeshElement* > & boundaryVolumes );
1069     void addSegments( SMESH_MesherHelper&                      helper,
1070                       const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap );
1071     void getVolumes( vector< const SMDS_MeshElement* > & volumes );
1072     void getBoundaryElems( vector< const SMDS_MeshElement* > & boundaryVolumes );
1073     void removeExcessSideDivision(const vector< Hexahedron* >& allHexa);
1074     void removeExcessNodes(vector< Hexahedron* >& allHexa);
1075     void preventVolumesOverlapping();
1076     TGeomID getAnyFace() const;
1077     void cutByExtendedInternal( std::vector< Hexahedron* >& hexes,
1078                                 const TColStd_MapOfInteger& intEdgeIDs );
1079     gp_Pnt mostDistantInternalPnt( int hexIndex, const gp_Pnt& p1, const gp_Pnt& p2 );
1080     bool isOutPoint( _Link& link, int iP, SMESH_MesherHelper& helper, const Solid* solid ) const;
1081     void sortVertexNodes(vector<_Node*>& nodes, _Node* curNode, TGeomID face);
1082     bool isInHole() const;
1083     bool hasStrangeEdge() const;
1084     bool checkPolyhedronSize( bool isCutByInternalFace, double & volSize ) const;
1085     int checkPolyhedronValidity( _volumeDef* volDef, std::vector<std::vector<int>>& splitQuantities, 
1086                                   std::vector<std::vector<const SMDS_MeshNode*>>& splitNodes );
1087     const SMDS_MeshElement* addPolyhedronToMesh( _volumeDef* volDef,  SMESH_MesherHelper& helper, const std::vector<const SMDS_MeshNode*>& nodes, 
1088                                                 const std::vector<int>& quantities );
1089     bool addHexa ();
1090     bool addTetra();
1091     bool addPenta();
1092     bool addPyra ();
1093     bool debugDumpLink( _Link* link );
1094     _Node* findEqualNode( vector< _Node* >&       nodes,
1095                           const E_IntersectPoint* ip,
1096                           const double            tol2 )
1097     {
1098       for ( size_t i = 0; i < nodes.size(); ++i )
1099         if ( nodes[i]->EdgeIntPnt() == ip ||
1100              nodes[i]->Point().SquareDistance( ip->_point ) <= tol2 )
1101           return nodes[i];
1102       return 0;
1103     }
1104     bool isCorner( const _Node* node ) const { return ( node >= &_hexNodes[0] &&
1105                                                         node -  &_hexNodes[0] < 8 ); }
1106     bool hasEdgesAround( const ConcaveFace* cf ) const;
1107     bool isImplementEdges() const { return _grid->_edgeIntPool.nbElements(); }
1108     bool isOutParam(const double uvw[3]) const;
1109
1110     typedef boost::container::flat_map< TGeomID, size_t > TID2Nb;
1111     static void insertAndIncrement( TGeomID id, TID2Nb& id2nbMap )
1112     {
1113       TID2Nb::value_type s0( id, 0 );
1114       TID2Nb::iterator id2nb = id2nbMap.insert( s0 ).first;
1115       id2nb->second++;
1116     }
1117   }; // class Hexahedron
1118
1119 #ifdef WITH_TBB
1120   // --------------------------------------------------------------------------
1121   /*!
1122    * \brief Hexahedron computing volumes in one thread
1123    */
1124   struct ParallelHexahedron
1125   {
1126     vector< Hexahedron* >& _hexVec;
1127     ParallelHexahedron( vector< Hexahedron* >& hv ): _hexVec(hv) {}
1128     void operator() ( const tbb::blocked_range<size_t>& r ) const
1129     {
1130       for ( size_t i = r.begin(); i != r.end(); ++i )
1131         if ( Hexahedron* hex = _hexVec[ i ] )
1132           hex->computeElements();
1133     }
1134   };
1135   // --------------------------------------------------------------------------
1136   /*!
1137    * \brief Structure intersecting certain nb of faces with GridLine's in one thread
1138    */
1139   struct ParallelIntersector
1140   {
1141     vector< FaceGridIntersector >& _faceVec;
1142     ParallelIntersector( vector< FaceGridIntersector >& faceVec): _faceVec(faceVec){}
1143     void operator() ( const tbb::blocked_range<size_t>& r ) const
1144     {
1145       for ( size_t i = r.begin(); i != r.end(); ++i )
1146         _faceVec[i].Intersect();
1147     }
1148   };
1149 #endif
1150
1151   //=============================================================================
1152   // Implementation of internal utils
1153   //=============================================================================
1154   /*!
1155    * \brief adjust \a i to have \a val between values[i] and values[i+1]
1156    */
1157   inline void locateValue( int & i, double val, const vector<double>& values,
1158                            int& di, double tol )
1159   {
1160     //val += values[0]; // input \a val is measured from 0.
1161     if ( i > (int) values.size()-2 )
1162       i = values.size()-2;
1163     else
1164       while ( i+2 < (int) values.size() && val > values[ i+1 ])
1165         ++i;
1166     while ( i > 0 && val < values[ i ])
1167       --i;
1168
1169     if ( i > 0 && val - values[ i ] < tol )
1170       di = -1;
1171     else if ( i+2 < (int) values.size() && values[ i+1 ] - val < tol )
1172       di = 1;
1173     else
1174       di = 0;
1175   }
1176   //=============================================================================
1177   /*
1178    * Return a vector of SOLIDS sharing given shapes
1179    */
1180   GeomIDVecHelder Geometry::GetSolidIDsByShapeID( const vector< TGeomID >& theShapeIDs ) const
1181   {
1182     if ( theShapeIDs.size() == 1 )
1183       return GeomIDVecHelder( _solidIDsByShapeID[ theShapeIDs[ 0 ]], /*owner=*/false );
1184
1185     // look for an empty slot in _solidIDsByShapeID
1186     vector< TGeomID > * resultIDs = 0;
1187     for ( const vector< TGeomID >& vec : _solidIDsByShapeID )
1188       if ( vec.empty() )
1189       {
1190         resultIDs = const_cast< vector< TGeomID > * >( & vec );
1191         break;
1192       }
1193     // fill in resultIDs
1194     for ( const TGeomID& id : theShapeIDs )
1195       for ( const TGeomID& solid : _solidIDsByShapeID[ id ])
1196       {
1197         if ( std::find( resultIDs->begin(), resultIDs->end(), solid ) == resultIDs->end() )
1198           resultIDs->push_back( solid );
1199       }
1200     return GeomIDVecHelder( *resultIDs, /*owner=*/true );
1201   }
1202   //=============================================================================
1203   /*
1204    * Remove coincident intersection points
1205    */
1206   void GridLine::RemoveExcessIntPoints( const double tol )
1207   {
1208     if ( _intPoints.size() < 2 ) return;
1209
1210     set< Transition > tranSet;
1211     multiset< F_IntersectPoint >::iterator ip1, ip2 = _intPoints.begin();
1212     while ( ip2 != _intPoints.end() )
1213     {
1214       tranSet.clear();
1215       ip1 = ip2++;
1216       while ( ip2 != _intPoints.end() && ip2->_paramOnLine - ip1->_paramOnLine <= tol )
1217       {
1218         tranSet.insert( ip1->_transition );
1219         tranSet.insert( ip2->_transition );
1220         ip2->Add( ip1->_faceIDs );
1221         _intPoints.erase( ip1 );
1222         ip1 = ip2++;
1223       }
1224       if ( tranSet.size() > 1 ) // points with different transition coincide
1225       {
1226         bool isIN  = tranSet.count( Trans_IN );
1227         bool isOUT = tranSet.count( Trans_OUT );
1228         if ( isIN && isOUT )
1229           (*ip1)._transition = Trans_TANGENT;
1230         else
1231           (*ip1)._transition = isIN ? Trans_IN : Trans_OUT;
1232       }
1233     }
1234   }
1235   //================================================================================
1236   /*
1237    * Return ID of SOLID for nodes before the given intersection point
1238    */
1239   TGeomID GridLine::GetSolidIDBefore( multiset< F_IntersectPoint >::iterator ip,
1240                                       const TGeomID                          prevID,
1241                                       const Geometry&                        geom )
1242   {
1243     if ( ip == _intPoints.begin() )
1244       return 0;
1245
1246     if ( geom.IsOneSolid() )
1247     {
1248       bool isOut = true;
1249       switch ( ip->_transition ) {
1250       case Trans_IN:      isOut = true;            break;
1251       case Trans_OUT:     isOut = false;           break;
1252       case Trans_TANGENT: isOut = ( prevID != 0 ); break;
1253       case Trans_APEX:
1254       {
1255         // singularity point (apex of a cone)
1256         multiset< F_IntersectPoint >::iterator ipBef = ip, ipAft = ++ip;
1257         if ( ipAft == _intPoints.end() )
1258           isOut = false;
1259         else
1260         {
1261           --ipBef;
1262           if ( ipBef->_transition != ipAft->_transition )
1263             isOut = ( ipBef->_transition == Trans_OUT );
1264           else
1265             isOut = ( ipBef->_transition != Trans_OUT );
1266         }
1267         break;
1268       }
1269       case Trans_INTERNAL: isOut = false;
1270       default:;
1271       }
1272       return isOut ? 0 : geom._soleSolid.ID();
1273     }
1274
1275     GeomIDVecHelder solids = geom.GetSolidIDsByShapeID( ip->_faceIDs );
1276
1277     --ip;
1278     if ( ip->_transition == Trans_INTERNAL )
1279       return prevID;
1280
1281     GeomIDVecHelder solidsBef = geom.GetSolidIDsByShapeID( ip->_faceIDs );
1282
1283     if ( ip->_transition == Trans_IN ||
1284          ip->_transition == Trans_OUT )
1285     {
1286       if ( solidsBef.size() == 1 )
1287       {
1288         if ( solidsBef[0] == prevID )
1289           return ip->_transition == Trans_OUT ? 0 : solidsBef[0];
1290         else
1291           return solidsBef[0];
1292       }
1293
1294       if ( solids.size() == 2 )
1295       {
1296         if ( solids == solidsBef )
1297           return solids.contain( prevID ) ? solids.otherThan( prevID ) : theUndefID; // bos #29212
1298       }
1299       return solids.oneCommon( solidsBef );
1300     }
1301
1302     if ( solidsBef.size() == 1 )
1303       return solidsBef[0];
1304
1305     return solids.oneCommon( solidsBef );
1306   }
1307   //================================================================================
1308   /*
1309    * Adds face IDs
1310    */
1311   bool B_IntersectPoint::Add( const vector< TGeomID >& fIDs,
1312                               const SMDS_MeshNode*     n) const
1313   {
1314     const std::lock_guard<std::mutex> lock(_bMutex);
1315     size_t prevNbF = _faceIDs.size();
1316
1317     if ( _faceIDs.empty() )
1318       _faceIDs = fIDs;
1319     else
1320       for ( size_t i = 0; i < fIDs.size(); ++i )
1321       {
1322         vector< TGeomID >::iterator it =
1323           std::find( _faceIDs.begin(), _faceIDs.end(), fIDs[i] );
1324         if ( it == _faceIDs.end() )
1325           _faceIDs.push_back( fIDs[i] );
1326       }
1327     if ( !_node && n != NULL )
1328       _node = n;
1329
1330     return prevNbF < _faceIDs.size();
1331   }
1332   //================================================================================
1333   /*
1334    * Return ID of a common face if any, else zero
1335    */
1336   TGeomID B_IntersectPoint::HasCommonFace( const B_IntersectPoint * other, TGeomID avoidFace ) const
1337   {
1338     if ( other )
1339       for ( size_t i = 0; i < other->_faceIDs.size(); ++i )
1340         if ( avoidFace != other->_faceIDs[i] &&
1341              IsOnFace   ( other->_faceIDs[i] ))
1342           return other->_faceIDs[i];
1343     return 0;
1344   }
1345   //================================================================================
1346   /*
1347    * Return faces common with other point
1348    */
1349   size_t B_IntersectPoint::GetCommonFaces( const B_IntersectPoint * other, TGeomID* common ) const
1350   {
1351     size_t nbComm = 0;
1352     if ( !other )
1353       return nbComm;
1354     if ( _faceIDs.size() > other->_faceIDs.size() )
1355       return other->GetCommonFaces( this, common );
1356     for ( const TGeomID& face : _faceIDs )
1357       if ( other->IsOnFace( face ))
1358         common[ nbComm++ ] = face;
1359     return nbComm;
1360   }
1361   //================================================================================
1362   /*
1363    * Return \c true if \a faceID in in this->_faceIDs
1364    */
1365   bool B_IntersectPoint::IsOnFace( TGeomID faceID ) const // returns true if faceID is found
1366   {
1367     vector< TGeomID >::const_iterator it =
1368       std::find( _faceIDs.begin(), _faceIDs.end(), faceID );
1369     return ( it != _faceIDs.end() );
1370   }
1371   //================================================================================
1372   /*
1373    * OneOfSolids initialization
1374    */
1375   void OneOfSolids::Init( const TopoDS_Shape& solid,
1376                           TopAbs_ShapeEnum    subType,
1377                           const SMESHDS_Mesh* mesh )
1378   {
1379     SetID( mesh->ShapeToIndex( solid ));
1380
1381     if ( subType == TopAbs_FACE )
1382       SetHasInternalFaces( false );
1383
1384     for ( TopExp_Explorer sub( solid, subType ); sub.More(); sub.Next() )
1385     {
1386       _subIDs.Add( mesh->ShapeToIndex( sub.Current() ));
1387       if ( subType == TopAbs_FACE )
1388       {
1389         _faces.Add( sub.Current() );
1390         if ( sub.Current().Orientation() == TopAbs_INTERNAL )
1391           SetHasInternalFaces( true );
1392
1393         TGeomID faceID = mesh->ShapeToIndex( sub.Current() );
1394         if ( sub.Current().Orientation() == TopAbs_INTERNAL ||
1395              sub.Current().Orientation() == mesh->IndexToShape( faceID ).Orientation() )
1396           _outFaceIDs.Add( faceID );
1397       }
1398     }
1399   }
1400   //================================================================================
1401   /*
1402    * Return an iterator on GridLine's in a given direction
1403    */
1404   LineIndexer Grid::GetLineIndexer(size_t iDir) const
1405   {
1406     const size_t indices[] = { 1,2,0, 0,2,1, 0,1,2 };
1407     const string s      [] = { "X", "Y", "Z" };
1408     LineIndexer li( _coords[0].size(),  _coords[1].size(),    _coords[2].size(),
1409                     indices[iDir*3],    indices[iDir*3+1],    indices[iDir*3+2],
1410                     s[indices[iDir*3]], s[indices[iDir*3+1]], s[indices[iDir*3+2]]);
1411     return li;
1412   }
1413   //================================================================================
1414   /*
1415    * Return direction [0,1,2] of a GridLine
1416    */
1417   size_t Grid::GetLineDir( const GridLine* line, size_t & index ) const
1418   {
1419     for ( size_t iDir = 0; iDir < 3; ++iDir )
1420       if ( &_lines[ iDir ][0] <= line && line <= &_lines[ iDir ].back() )
1421       {
1422         index = line - &_lines[ iDir ][0];
1423         return iDir;
1424       }
1425     return -1;
1426   }
1427   //=============================================================================
1428   /*
1429    * Creates GridLine's of the grid
1430    */
1431   void Grid::SetCoordinates(const vector<double>& xCoords,
1432                             const vector<double>& yCoords,
1433                             const vector<double>& zCoords,
1434                             const double*         axesDirs,
1435                             const Bnd_Box&        shapeBox)
1436   {
1437     _coords[0] = xCoords;
1438     _coords[1] = yCoords;
1439     _coords[2] = zCoords;
1440
1441     _axes[0].SetCoord( axesDirs[0],
1442                        axesDirs[1],
1443                        axesDirs[2]);
1444     _axes[1].SetCoord( axesDirs[3],
1445                        axesDirs[4],
1446                        axesDirs[5]);
1447     _axes[2].SetCoord( axesDirs[6],
1448                        axesDirs[7],
1449                        axesDirs[8]);
1450     _axes[0].Normalize();
1451     _axes[1].Normalize();
1452     _axes[2].Normalize();
1453
1454     _invB.SetCols( _axes[0], _axes[1], _axes[2] );
1455     _invB.Invert();
1456
1457     // compute tolerance
1458     _minCellSize = Precision::Infinite();
1459     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
1460     {
1461       for ( size_t i = 1; i < _coords[ iDir ].size(); ++i )
1462       {
1463         double cellLen = _coords[ iDir ][ i ] - _coords[ iDir ][ i-1 ];
1464         if ( cellLen < _minCellSize )
1465           _minCellSize = cellLen;
1466       }
1467     }
1468     if ( _minCellSize < Precision::Confusion() )
1469       throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
1470                                 SMESH_Comment("Too small cell size: ") << _minCellSize );
1471     _tol = _minCellSize / 1000.;
1472
1473     // attune grid extremities to shape bounding box
1474
1475     double sP[6]; // aXmin, aYmin, aZmin, aXmax, aYmax, aZmax
1476     shapeBox.Get(sP[0],sP[1],sP[2],sP[3],sP[4],sP[5]);
1477     double* cP[6] = { &_coords[0].front(), &_coords[1].front(), &_coords[2].front(),
1478                       &_coords[0].back(),  &_coords[1].back(),  &_coords[2].back() };
1479     for ( int i = 0; i < 6; ++i )
1480       if ( fabs( sP[i] - *cP[i] ) < _tol )
1481         *cP[i] = sP[i];// + _tol/1000. * ( i < 3 ? +1 : -1 );
1482
1483     for ( int iDir = 0; iDir < 3; ++iDir )
1484     {
1485       if ( _coords[iDir][0] - sP[iDir] > _tol )
1486       {
1487         _minCellSize = Min( _minCellSize, _coords[iDir][0] - sP[iDir] );
1488         _coords[iDir].insert( _coords[iDir].begin(), sP[iDir] + _tol/1000.);
1489       }
1490       if ( sP[iDir+3] - _coords[iDir].back() > _tol  )
1491       {
1492         _minCellSize = Min( _minCellSize, sP[iDir+3] - _coords[iDir].back() );
1493         _coords[iDir].push_back( sP[iDir+3] - _tol/1000.);
1494       }
1495     }
1496     _tol = _minCellSize / 1000.;
1497
1498     _origin = ( _coords[0][0] * _axes[0] +
1499                 _coords[1][0] * _axes[1] +
1500                 _coords[2][0] * _axes[2] );
1501
1502     // create lines
1503     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
1504     {
1505       LineIndexer li = GetLineIndexer( iDir );
1506       _lines[iDir].resize( li.NbLines() );
1507       double len = _coords[ iDir ].back() - _coords[iDir].front();
1508       for ( ; li.More(); ++li )
1509       {
1510         GridLine& gl = _lines[iDir][ li.LineIndex() ];
1511         gl._line.SetLocation( _coords[0][li.I()] * _axes[0] +
1512                               _coords[1][li.J()] * _axes[1] +
1513                               _coords[2][li.K()] * _axes[2] );
1514         gl._line.SetDirection( _axes[ iDir ]);
1515         gl._length = len;
1516       }
1517     }
1518   }
1519   //================================================================================
1520   /*
1521    * Return local ID of shape
1522    */
1523   TGeomID Grid::ShapeID( const TopoDS_Shape& s ) const
1524   {
1525     return _helper->GetMeshDS()->ShapeToIndex( s );
1526   }
1527   //================================================================================
1528   /*
1529    * Return a shape by its local ID
1530    */
1531   const TopoDS_Shape& Grid::Shape( TGeomID id ) const
1532   {
1533     return _helper->GetMeshDS()->IndexToShape( id );
1534   }
1535   //================================================================================
1536   /*
1537    * Initialize _geometry
1538    */
1539   void Grid::InitGeometry( const TopoDS_Shape& theShapeToMesh )
1540   {
1541     SMESH_Mesh* mesh = _helper->GetMesh();
1542
1543     _geometry._mainShape = theShapeToMesh;
1544     _geometry._extIntFaceID = mesh->GetMeshDS()->MaxShapeIndex() * 100;
1545     _geometry._soleSolid.SetID( 0 );
1546     _geometry._soleSolid.SetHasInternalFaces( false );
1547
1548     InitClassifier( theShapeToMesh, TopAbs_VERTEX, _geometry._vertexClassifier );
1549     InitClassifier( theShapeToMesh, TopAbs_EDGE  , _geometry._edgeClassifier );
1550
1551     TopExp_Explorer solidExp( theShapeToMesh, TopAbs_SOLID );
1552
1553     bool isSeveralSolids = false;
1554     if ( _toConsiderInternalFaces ) // check nb SOLIDs
1555     {
1556       solidExp.Next();
1557       isSeveralSolids = solidExp.More();
1558       _toConsiderInternalFaces = isSeveralSolids;
1559       solidExp.ReInit();
1560
1561       if ( !isSeveralSolids ) // look for an internal FACE
1562       {
1563         TopExp_Explorer fExp( theShapeToMesh, TopAbs_FACE );
1564         for ( ; fExp.More() &&  !_toConsiderInternalFaces; fExp.Next() )
1565           _toConsiderInternalFaces = ( fExp.Current().Orientation() == TopAbs_INTERNAL );
1566
1567         _geometry._soleSolid.SetHasInternalFaces( _toConsiderInternalFaces );
1568         _geometry._soleSolid.SetID( ShapeID( solidExp.Current() ));
1569       }
1570       else // fill Geometry::_solidByID
1571       {
1572         for ( ; solidExp.More(); solidExp.Next() )
1573         {
1574           OneOfSolids & solid = _geometry._solidByID[ ShapeID( solidExp.Current() )];
1575           solid.Init( solidExp.Current(), TopAbs_FACE,   mesh->GetMeshDS() );
1576           solid.Init( solidExp.Current(), TopAbs_EDGE,   mesh->GetMeshDS() );
1577           solid.Init( solidExp.Current(), TopAbs_VERTEX, mesh->GetMeshDS() );
1578         }
1579       }
1580     }
1581     else
1582     {
1583       _geometry._soleSolid.SetID( ShapeID( solidExp.Current() ));
1584     }
1585
1586     if ( !_toCreateFaces )
1587     {
1588       int nbSolidsGlobal = _helper->Count( mesh->GetShapeToMesh(), TopAbs_SOLID, false );
1589       int nbSolidsLocal  = _helper->Count( theShapeToMesh,         TopAbs_SOLID, false );
1590       _toCreateFaces = ( nbSolidsLocal < nbSolidsGlobal );
1591     }
1592
1593     TopTools_IndexedMapOfShape faces;
1594     TopExp::MapShapes( theShapeToMesh, TopAbs_FACE, faces );
1595
1596     // find boundary FACEs on boundary of mesh->ShapeToMesh()
1597     if ( _toCreateFaces )
1598       for ( int i = 1; i <= faces.Size(); ++i )
1599         if ( faces(i).Orientation() != TopAbs_INTERNAL &&
1600              _helper->NbAncestors( faces(i), *mesh, TopAbs_SOLID ) == 1 )
1601         {
1602           _geometry._boundaryFaces.Add( ShapeID( faces(i) ));
1603         }
1604
1605     if ( isSeveralSolids )
1606       for ( int i = 1; i <= faces.Size(); ++i )
1607       {
1608         SetSolidFather( faces(i), theShapeToMesh );
1609         for ( TopExp_Explorer eExp( faces(i), TopAbs_EDGE ); eExp.More(); eExp.Next() )
1610         {
1611           const TopoDS_Edge& edge = TopoDS::Edge( eExp.Current() );
1612           SetSolidFather( edge, theShapeToMesh );
1613           SetSolidFather( _helper->IthVertex( 0, edge ), theShapeToMesh );
1614           SetSolidFather( _helper->IthVertex( 1, edge ), theShapeToMesh );
1615         }
1616       }
1617
1618     // fill in _geometry._shape2NbNodes == find already meshed sub-shapes
1619     _geometry._shape2NbNodes.Clear();
1620     if ( mesh->NbNodes() > 0 )
1621     {
1622       for ( TopAbs_ShapeEnum type : { TopAbs_FACE, TopAbs_EDGE, TopAbs_VERTEX })
1623         for ( TopExp_Explorer exp( theShapeToMesh, type ); exp.More(); exp.Next() )
1624         {
1625           if ( _geometry._shape2NbNodes.IsBound( exp.Current() ))
1626             continue;
1627           if ( SMESHDS_SubMesh* sm = mesh->GetMeshDS()->MeshElements( exp.Current() ))
1628             if ( sm->NbNodes() > 0 )
1629               _geometry._shape2NbNodes.Bind( exp.Current(), sm->NbNodes() );
1630         }
1631     }
1632
1633     // fill in Solid::_concaveVertex
1634     vector< TGeomID > soleSolidID( 1, _geometry._soleSolid.ID() );
1635     for ( int i = 1; i <= faces.Size(); ++i )
1636     {
1637       const TopoDS_Face& F = TopoDS::Face( faces( i ));
1638       TError error;
1639       TSideVector wires = StdMeshers_FaceSide::GetFaceWires( F, *mesh, 0, error,
1640                                                              nullptr, nullptr, false );
1641       for ( StdMeshers_FaceSidePtr& wire : wires )
1642       {
1643         const int nbEdges = wire->NbEdges();
1644         if ( nbEdges < 2 && SMESH_Algo::isDegenerated( wire->Edge(0)))
1645           continue;
1646         for ( int iE1 = 0; iE1 < nbEdges; ++iE1 )
1647         {
1648           if ( SMESH_Algo::isDegenerated( wire->Edge( iE1 ))) continue;
1649           int iE2 = ( iE1 + 1 ) % nbEdges;
1650           while ( SMESH_Algo::isDegenerated( wire->Edge( iE2 )))
1651             iE2 = ( iE2 + 1 ) % nbEdges;
1652           TopoDS_Vertex V = wire->FirstVertex( iE2 );
1653           double angle = _helper->GetAngle( wire->Edge( iE1 ),
1654                                             wire->Edge( iE2 ), F, V );
1655           if ( angle < -5. * M_PI / 180. )
1656           {
1657             TGeomID faceID = ShapeID( F );
1658             const vector< TGeomID > & solids =
1659               _geometry.IsOneSolid() ? soleSolidID : GetSolidIDs( faceID );
1660             for ( const TGeomID & solidID : solids )
1661             {
1662               Solid* solid = GetSolid( solidID );
1663               TGeomID   V1 = ShapeID( wire->FirstVertex( iE1 ));
1664               TGeomID   V2 = ShapeID( wire->LastVertex ( iE2 ));
1665               solid->SetConcave( ShapeID( V ), faceID,
1666                                  wire->EdgeID( iE1 ), wire->EdgeID( iE2 ), V1, V2 );
1667             }
1668           }
1669         }
1670       }
1671     }
1672
1673     return;
1674   }
1675   //================================================================================
1676   /*
1677    * Store ID of SOLID as father of its child shape ID
1678    */
1679   void Grid::SetSolidFather( const TopoDS_Shape& s, const TopoDS_Shape& theShapeToMesh )
1680   {
1681     if ( _geometry._solidIDsByShapeID.empty() )
1682       _geometry._solidIDsByShapeID.resize( _helper->GetMeshDS()->MaxShapeIndex() + 1 );
1683
1684     vector< TGeomID > & solidIDs = _geometry._solidIDsByShapeID[ ShapeID( s )];
1685     if ( !solidIDs.empty() )
1686       return;
1687     solidIDs.reserve(2);
1688     PShapeIteratorPtr solidIt = _helper->GetAncestors( s,
1689                                                        *_helper->GetMesh(),
1690                                                        TopAbs_SOLID,
1691                                                        & theShapeToMesh );
1692     while ( const TopoDS_Shape* solid = solidIt->next() )
1693       solidIDs.push_back( ShapeID( *solid ));
1694   }
1695   //================================================================================
1696   /*
1697    * Return IDs of solids given sub-shape belongs to
1698    */
1699   const vector< TGeomID > & Grid::GetSolidIDs( TGeomID subShapeID ) const
1700   {
1701     return _geometry._solidIDsByShapeID[ subShapeID ];
1702   }
1703   //================================================================================
1704   /*
1705    * Check if a sub-shape belongs to several SOLIDs
1706    */
1707   bool Grid::IsShared( TGeomID shapeID ) const
1708   {
1709     return !_geometry.IsOneSolid() && ( _geometry._solidIDsByShapeID[ shapeID ].size() > 1 );
1710   }
1711   //================================================================================
1712   /*
1713    * Check if any of FACEs belongs to several SOLIDs
1714    */
1715   bool Grid::IsAnyShared( const std::vector< TGeomID >& faceIDs ) const
1716   {
1717     for ( size_t i = 0; i < faceIDs.size(); ++i )
1718       if ( IsShared( faceIDs[ i ]))
1719         return true;
1720     return false;
1721   }
1722   //================================================================================
1723   /*
1724    * Return Solid by ID
1725    */
1726   Solid* Grid::GetSolid( TGeomID solidID )
1727   {
1728     if ( !solidID || _geometry.IsOneSolid() || _geometry._solidByID.empty() )
1729       return & _geometry._soleSolid;
1730
1731     return & _geometry._solidByID[ solidID ];
1732   }
1733   //================================================================================
1734   /*
1735    * Return OneOfSolids by ID
1736    */
1737   Solid* Grid::GetOneOfSolids( TGeomID solidID )
1738   {
1739     map< TGeomID, OneOfSolids >::iterator is2s = _geometry._solidByID.find( solidID );
1740     if ( is2s != _geometry._solidByID.end() )
1741       return & is2s->second;
1742
1743     return & _geometry._soleSolid;
1744   }
1745   //================================================================================
1746   /*
1747    * Check if transition on given FACE is correct for a given SOLID
1748    */
1749   bool Grid::IsCorrectTransition( TGeomID faceID, const Solid* solid )
1750   {
1751     if ( _geometry.IsOneSolid() )
1752       return true;
1753
1754     const vector< TGeomID >& solidIDs = _geometry._solidIDsByShapeID[ faceID ];
1755     return solidIDs[0] == solid->ID();
1756   }
1757
1758   //================================================================================
1759   /*
1760    * Assign to geometry a node at FACE intersection
1761    * Return a found supporting VERTEX
1762    */
1763   void Grid::SetOnShape( const SMDS_MeshNode* n, const F_IntersectPoint& ip,
1764                          TopoDS_Vertex* vertex, bool unset )
1765   {
1766     TopoDS_Shape s;
1767     SMESHDS_Mesh* mesh = _helper->GetMeshDS();
1768     if ( ip._faceIDs.size() == 1 )
1769     {
1770       mesh->SetNodeOnFace( n, ip._faceIDs[0], ip._u, ip._v );
1771     }
1772     else if ( _geometry._vertexClassifier.IsSatisfy( n, &s ))
1773     {
1774       if ( unset ) mesh->UnSetNodeOnShape( n );
1775       mesh->SetNodeOnVertex( n, TopoDS::Vertex( s ));
1776       if ( vertex )
1777         *vertex = TopoDS::Vertex( s );
1778     }
1779     else if ( _geometry._edgeClassifier.IsSatisfy( n, &s ))
1780     {
1781       if ( unset ) mesh->UnSetNodeOnShape( n );
1782       mesh->SetNodeOnEdge( n, TopoDS::Edge( s ));
1783     }
1784     else if ( ip._faceIDs.size() > 0 )
1785     {
1786       mesh->SetNodeOnFace( n, ip._faceIDs[0], ip._u, ip._v );
1787     }
1788     else if ( !unset && _geometry.IsOneSolid() )
1789     {
1790       mesh->SetNodeInVolume( n, _geometry._soleSolid.ID() );
1791     }
1792   }
1793   //================================================================================
1794   /*
1795    * Fill in B_IntersectPoint::_faceIDs with all FACEs sharing a VERTEX
1796    */
1797   void Grid::UpdateFacesOfVertex( const B_IntersectPoint& ip, const TopoDS_Vertex& vertex )
1798   {
1799     if ( vertex.IsNull() )
1800       return;
1801     std::vector< int > faceID(1);
1802     PShapeIteratorPtr fIt = _helper->GetAncestors( vertex, *_helper->GetMesh(),
1803                                                    TopAbs_FACE, & _geometry._mainShape );
1804     while ( const TopoDS_Shape* face = fIt->next() )
1805     {
1806       faceID[ 0 ] = ShapeID( *face );
1807       ip.Add( faceID );
1808     }
1809   }
1810   //================================================================================
1811   /*
1812    * Initialize a classifier
1813    */
1814   void Grid::InitClassifier( const TopoDS_Shape&        mainShape,
1815                              TopAbs_ShapeEnum           shapeType,
1816                              Controls::ElementsOnShape& classifier )
1817   {
1818     TopTools_IndexedMapOfShape shapes;
1819     TopExp::MapShapes( mainShape, shapeType, shapes );
1820
1821     TopoDS_Compound compound; BRep_Builder builder;
1822     builder.MakeCompound( compound );
1823     for ( int i = 1; i <= shapes.Size(); ++i )
1824       builder.Add( compound, shapes(i) );
1825
1826     classifier.SetMesh( _helper->GetMeshDS() );
1827     //classifier.SetTolerance( _tol ); // _tol is not initialised
1828     classifier.SetShape( compound, SMDSAbs_Node );
1829   }
1830
1831   //================================================================================
1832   /*
1833    * Return EDGEs with FACEs to implement into the mesh
1834    */
1835   void Grid::GetEdgesToImplement( map< TGeomID, vector< TGeomID > > & edge2faceIDsMap,
1836                                   const TopoDS_Shape&                 shape,
1837                                   const vector< TopoDS_Shape >&       faces )
1838   {
1839     // check if there are strange EDGEs
1840     TopTools_IndexedMapOfShape faceMap;
1841     TopExp::MapShapes( _helper->GetMesh()->GetShapeToMesh(), TopAbs_FACE, faceMap );
1842     int nbFacesGlobal = faceMap.Size();
1843     faceMap.Clear( false );
1844     TopExp::MapShapes( shape, TopAbs_FACE, faceMap );
1845     int nbFacesLocal  = faceMap.Size();
1846     bool hasStrangeEdges = ( nbFacesGlobal > nbFacesLocal );
1847     if ( !_toAddEdges && !hasStrangeEdges )
1848       return; // no FACEs in contact with those meshed by other algo
1849
1850     for ( size_t i = 0; i < faces.size(); ++i )
1851     {
1852       _helper->SetSubShape( faces[i] );
1853       for ( TopExp_Explorer eExp( faces[i], TopAbs_EDGE ); eExp.More(); eExp.Next() )
1854       {
1855         const TopoDS_Edge& edge = TopoDS::Edge( eExp.Current() );
1856         if ( hasStrangeEdges )
1857         {
1858           bool hasStrangeFace = false;
1859           PShapeIteratorPtr faceIt = _helper->GetAncestors( edge, *_helper->GetMesh(), TopAbs_FACE);
1860           while ( const TopoDS_Shape* face = faceIt->next() )
1861             if (( hasStrangeFace = !faceMap.Contains( *face )))
1862               break;
1863           if ( !hasStrangeFace && !_toAddEdges )
1864             continue;
1865           _geometry._strangeEdges.Add( ShapeID( edge ));
1866           _geometry._strangeEdges.Add( ShapeID( _helper->IthVertex( 0, edge )));
1867           _geometry._strangeEdges.Add( ShapeID( _helper->IthVertex( 1, edge )));
1868         }
1869         if ( !SMESH_Algo::isDegenerated( edge ) &&
1870              !_helper->IsRealSeam( edge ))
1871         {
1872           edge2faceIDsMap[ ShapeID( edge )].push_back( ShapeID( faces[i] ));
1873         }
1874       }
1875     }
1876     return;
1877   }
1878
1879   //================================================================================
1880   /*
1881    * Computes coordinates of a point in the grid CS
1882    */
1883   void Grid::ComputeUVW(const gp_XYZ& P, double UVW[3])
1884   {
1885     gp_XYZ p = P * _invB;
1886     p.Coord( UVW[0], UVW[1], UVW[2] );
1887   }
1888   //================================================================================
1889   /*
1890    * Creates all nodes
1891    */
1892   void Grid::ComputeNodes(SMESH_MesherHelper& helper)
1893   {
1894     // state of each node of the grid relative to the geometry
1895     const size_t nbGridNodes = _coords[0].size() * _coords[1].size() * _coords[2].size();
1896     vector< TGeomID > shapeIDVec( nbGridNodes, theUndefID );
1897     _nodes.resize( nbGridNodes, 0 );
1898     _allBorderNodes.resize( nbGridNodes, 0 );
1899     _gridIntP.resize( nbGridNodes, NULL );
1900
1901     SMESHDS_Mesh* mesh = helper.GetMeshDS();
1902
1903     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
1904     {
1905       LineIndexer li = GetLineIndexer( iDir );
1906
1907       // find out a shift of node index while walking along a GridLine in this direction
1908       li.SetIndexOnLine( 0 );
1909       size_t nIndex0 = NodeIndex( li.I(), li.J(), li.K() );
1910       li.SetIndexOnLine( 1 );
1911       const size_t nShift = NodeIndex( li.I(), li.J(), li.K() ) - nIndex0;
1912       
1913       const vector<double> & coords = _coords[ iDir ];
1914       for ( ; li.More(); ++li ) // loop on lines in iDir
1915       {
1916         li.SetIndexOnLine( 0 );
1917         nIndex0 = NodeIndex( li.I(), li.J(), li.K() );
1918
1919         GridLine& line = _lines[ iDir ][ li.LineIndex() ];
1920         const gp_XYZ lineLoc = line._line.Location().XYZ();
1921         const gp_XYZ lineDir = line._line.Direction().XYZ();
1922
1923         line.RemoveExcessIntPoints( _tol );
1924         multiset< F_IntersectPoint >&     intPnts = line._intPoints;
1925         multiset< F_IntersectPoint >::iterator ip = intPnts.begin();
1926
1927         // Create mesh nodes at intersections with geometry
1928         // and set OUT state of nodes between intersections
1929
1930         TGeomID solidID = 0;
1931         const double* nodeCoord = & coords[0];
1932         const double* coord0    = nodeCoord;
1933         const double* coordEnd  = coord0 + coords.size();
1934         double nodeParam = 0;
1935         for ( ; ip != intPnts.end(); ++ip )
1936         {
1937           solidID = line.GetSolidIDBefore( ip, solidID, _geometry );
1938
1939           // set OUT state or just skip IN nodes before ip
1940           if ( nodeParam < ip->_paramOnLine - _tol )
1941           {
1942             while ( nodeParam < ip->_paramOnLine - _tol )
1943             {
1944               TGeomID & nodeShapeID = shapeIDVec[ nIndex0 + nShift * ( nodeCoord-coord0 ) ];
1945               nodeShapeID = Min( solidID, nodeShapeID );
1946               if ( ++nodeCoord <  coordEnd )
1947                 nodeParam = *nodeCoord - *coord0;
1948               else
1949                 break;                
1950             }
1951             if ( nodeCoord == coordEnd ) break;
1952           }
1953           
1954           // create a mesh node on a GridLine at ip if it does not coincide with a grid node
1955           if ( nodeParam > ip->_paramOnLine + _tol )
1956           {
1957             gp_XYZ xyz = lineLoc + ip->_paramOnLine * lineDir;
1958             ip->_node = mesh->AddNode( xyz.X(), xyz.Y(), xyz.Z() );
1959             ip->_indexOnLine = nodeCoord-coord0-1;
1960             TopoDS_Vertex v;
1961             SetOnShape( ip->_node, *ip, & v );
1962             UpdateFacesOfVertex( *ip, v );
1963           }
1964           // create a mesh node at ip coincident with a grid node
1965           else
1966           {
1967             int nodeIndex = nIndex0 + nShift * ( nodeCoord-coord0 );
1968             if ( !_nodes[ nodeIndex ] )
1969             {
1970               gp_XYZ xyz = lineLoc + nodeParam * lineDir;
1971               _nodes   [ nodeIndex ] = mesh->AddNode( xyz.X(), xyz.Y(), xyz.Z() );
1972               //_gridIntP[ nodeIndex ] = & * ip;
1973               //SetOnShape( _nodes[ nodeIndex ], *ip );
1974             }
1975             if ( _gridIntP[ nodeIndex ] )
1976               _gridIntP[ nodeIndex ]->Add( ip->_faceIDs );
1977             else
1978               _gridIntP[ nodeIndex ] = & * ip;
1979             // ip->_node        = _nodes[ nodeIndex ]; -- to differ from ip on links
1980             ip->_indexOnLine = nodeCoord-coord0;
1981             if ( ++nodeCoord < coordEnd )
1982               nodeParam = *nodeCoord - *coord0;
1983           }
1984         }
1985         // set OUT state to nodes after the last ip
1986         for ( ; nodeCoord < coordEnd; ++nodeCoord )
1987           shapeIDVec[ nIndex0 + nShift * ( nodeCoord-coord0 ) ] = 0;
1988       }
1989     }
1990
1991     // Create mesh nodes at !OUT nodes of the grid
1992
1993     for ( size_t z = 0; z < _coords[2].size(); ++z )
1994       for ( size_t y = 0; y < _coords[1].size(); ++y )
1995         for ( size_t x = 0; x < _coords[0].size(); ++x )
1996         {
1997           size_t nodeIndex = NodeIndex( x, y, z );
1998           if ( !_nodes[ nodeIndex ] &&
1999                0 < shapeIDVec[ nodeIndex ] && shapeIDVec[ nodeIndex ] < theUndefID )
2000           {
2001             gp_XYZ xyz = ( _coords[0][x] * _axes[0] +
2002                            _coords[1][y] * _axes[1] +
2003                            _coords[2][z] * _axes[2] );
2004             _nodes[ nodeIndex ] = mesh->AddNode( xyz.X(), xyz.Y(), xyz.Z() );
2005             mesh->SetNodeInVolume( _nodes[ nodeIndex ], shapeIDVec[ nodeIndex ]);
2006           }
2007           else if ( _nodes[ nodeIndex ] && _gridIntP[ nodeIndex ] /*&&
2008                     !_nodes[ nodeIndex]->GetShapeID()*/ )
2009           {
2010             TopoDS_Vertex v;
2011             SetOnShape( _nodes[ nodeIndex ], *_gridIntP[ nodeIndex ], & v );
2012             UpdateFacesOfVertex( *_gridIntP[ nodeIndex ], v );
2013           }
2014           else if ( _toUseQuanta && !_allBorderNodes[ nodeIndex ] /*add all nodes outside the body. Used to reconstruct the hexahedrals when polys are not desired!*/)
2015           {
2016             gp_XYZ xyz = ( _coords[0][x] * _axes[0] +
2017                            _coords[1][y] * _axes[1] +
2018                            _coords[2][z] * _axes[2] );
2019             _allBorderNodes[ nodeIndex ] = mesh->AddNode( xyz.X(), xyz.Y(), xyz.Z() );
2020             mesh->SetNodeInVolume( _allBorderNodes[ nodeIndex ], shapeIDVec[ nodeIndex ]);
2021           }
2022         }
2023
2024 #ifdef _MY_DEBUG_
2025     // check validity of transitions
2026     const char* trName[] = { "TANGENT", "IN", "OUT", "APEX" };
2027     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
2028     {
2029       LineIndexer li = GetLineIndexer( iDir );
2030       for ( ; li.More(); ++li )
2031       {
2032         multiset< F_IntersectPoint >& intPnts = _lines[ iDir ][ li.LineIndex() ]._intPoints;
2033         if ( intPnts.empty() ) continue;
2034         if ( intPnts.size() == 1 )
2035         {
2036           if ( intPnts.begin()->_transition != Trans_TANGENT &&
2037                intPnts.begin()->_transition != Trans_APEX )
2038             throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
2039                                       SMESH_Comment("Wrong SOLE transition of GridLine (")
2040                                       << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
2041                                       << ") along " << li._nameConst
2042                                       << ": " << trName[ intPnts.begin()->_transition] );
2043         }
2044         else
2045         {
2046           if ( intPnts.begin()->_transition == Trans_OUT )
2047             throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
2048                                       SMESH_Comment("Wrong START transition of GridLine (")
2049                                       << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
2050                                       << ") along " << li._nameConst
2051                                       << ": " << trName[ intPnts.begin()->_transition ]);
2052           if ( intPnts.rbegin()->_transition == Trans_IN )
2053             throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
2054                                       SMESH_Comment("Wrong END transition of GridLine (")
2055                                       << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
2056                                       << ") along " << li._nameConst
2057                                       << ": " << trName[ intPnts.rbegin()->_transition ]);
2058         }
2059       }
2060     }
2061 #endif
2062
2063     return;
2064   }
2065
2066   //=============================================================================
2067   /*
2068    * Intersects TopoDS_Face with all GridLine's
2069    */
2070   void FaceGridIntersector::Intersect()
2071   {
2072     FaceLineIntersector intersector;
2073     intersector._surfaceInt = GetCurveFaceIntersector();
2074     intersector._tol        = _grid->_tol;
2075     intersector._transOut   = _face.Orientation() == TopAbs_REVERSED ? Trans_IN : Trans_OUT;
2076     intersector._transIn    = _face.Orientation() == TopAbs_REVERSED ? Trans_OUT : Trans_IN;
2077
2078     typedef void (FaceLineIntersector::* PIntFun )(const GridLine& gridLine);
2079     PIntFun interFunction;
2080
2081     bool isDirect = true;
2082     BRepAdaptor_Surface surf( _face );
2083     switch ( surf.GetType() ) {
2084     case GeomAbs_Plane:
2085       intersector._plane = surf.Plane();
2086       interFunction = &FaceLineIntersector::IntersectWithPlane;
2087       isDirect = intersector._plane.Direct();
2088       break;
2089     case GeomAbs_Cylinder:
2090       intersector._cylinder = surf.Cylinder();
2091       interFunction = &FaceLineIntersector::IntersectWithCylinder;
2092       isDirect = intersector._cylinder.Direct();
2093       break;
2094     case GeomAbs_Cone:
2095       intersector._cone = surf.Cone();
2096       interFunction = &FaceLineIntersector::IntersectWithCone;
2097       //isDirect = intersector._cone.Direct();
2098       break;
2099     case GeomAbs_Sphere:
2100       intersector._sphere = surf.Sphere();
2101       interFunction = &FaceLineIntersector::IntersectWithSphere;
2102       isDirect = intersector._sphere.Direct();
2103       break;
2104     case GeomAbs_Torus:
2105       intersector._torus = surf.Torus();
2106       interFunction = &FaceLineIntersector::IntersectWithTorus;
2107       //isDirect = intersector._torus.Direct();
2108       break;
2109     default:
2110       interFunction = &FaceLineIntersector::IntersectWithSurface;
2111     }
2112     if ( !isDirect )
2113       std::swap( intersector._transOut, intersector._transIn );
2114
2115     _intersections.clear();
2116     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
2117     {
2118       if ( surf.GetType() == GeomAbs_Plane )
2119       {
2120         // check if all lines in this direction are parallel to a plane
2121         if ( intersector._plane.Axis().IsNormal( _grid->_lines[iDir][0]._line.Position(),
2122                                                  Precision::Angular()))
2123           continue;
2124         // find out a transition, that is the same for all lines of a direction
2125         gp_Dir plnNorm = intersector._plane.Axis().Direction();
2126         gp_Dir lineDir = _grid->_lines[iDir][0]._line.Direction();
2127         intersector._transition =
2128           ( plnNorm * lineDir < 0 ) ? intersector._transIn : intersector._transOut;
2129       }
2130       if ( surf.GetType() == GeomAbs_Cylinder )
2131       {
2132         // check if all lines in this direction are parallel to a cylinder
2133         if ( intersector._cylinder.Axis().IsParallel( _grid->_lines[iDir][0]._line.Position(),
2134                                                       Precision::Angular()))
2135           continue;
2136       }
2137
2138       // intersect the grid lines with the face
2139       for ( size_t iL = 0; iL < _grid->_lines[iDir].size(); ++iL )
2140       {
2141         GridLine& gridLine = _grid->_lines[iDir][iL];
2142         if ( _bndBox.IsOut( gridLine._line )) continue;
2143
2144         intersector._intPoints.clear();
2145         (intersector.*interFunction)( gridLine ); // <- intersection with gridLine
2146         for ( size_t i = 0; i < intersector._intPoints.size(); ++i )
2147           _intersections.push_back( make_pair( &gridLine, intersector._intPoints[i] ));
2148       }
2149     }
2150
2151     if ( _face.Orientation() == TopAbs_INTERNAL )
2152     {
2153       for ( size_t i = 0; i < _intersections.size(); ++i )
2154         if ( _intersections[i].second._transition == Trans_IN ||
2155              _intersections[i].second._transition == Trans_OUT )
2156         {
2157           _intersections[i].second._transition = Trans_INTERNAL;
2158         }
2159     }
2160     return;
2161   }
2162   //================================================================================
2163   /*
2164    * Return true if (_u,_v) is on the face
2165    */
2166   bool FaceLineIntersector::UVIsOnFace() const
2167   {
2168     TopAbs_State state = _surfaceInt->ClassifyUVPoint(gp_Pnt2d( _u,_v ));
2169     return ( state == TopAbs_IN || state == TopAbs_ON );
2170   }
2171   //================================================================================
2172   /*
2173    * Store an intersection if it is IN or ON the face
2174    */
2175   void FaceLineIntersector::addIntPoint(const bool toClassify)
2176   {
2177     if ( !toClassify || UVIsOnFace() )
2178     {
2179       F_IntersectPoint p;
2180       p._paramOnLine = _w;
2181       p._u           = _u;
2182       p._v           = _v;
2183       p._transition  = _transition;
2184       _intPoints.push_back( p );
2185     }
2186   }
2187   //================================================================================
2188   /*
2189    * Intersect a line with a plane
2190    */
2191   void FaceLineIntersector::IntersectWithPlane(const GridLine& gridLine)
2192   {
2193     IntAna_IntConicQuad linPlane( gridLine._line, _plane, Precision::Angular());
2194     _w = linPlane.ParamOnConic(1);
2195     if ( isParamOnLineOK( gridLine._length ))
2196     {
2197       ElSLib::Parameters(_plane, linPlane.Point(1) ,_u,_v);
2198       addIntPoint();
2199     }
2200   }
2201   //================================================================================
2202   /*
2203    * Intersect a line with a cylinder
2204    */
2205   void FaceLineIntersector::IntersectWithCylinder(const GridLine& gridLine)
2206   {
2207     IntAna_IntConicQuad linCylinder( gridLine._line, _cylinder );
2208     if ( linCylinder.IsDone() && linCylinder.NbPoints() > 0 )
2209     {
2210       _w = linCylinder.ParamOnConic(1);
2211       if ( linCylinder.NbPoints() == 1 )
2212         _transition = Trans_TANGENT;
2213       else
2214         _transition = _w < linCylinder.ParamOnConic(2) ? _transIn : _transOut;
2215       if ( isParamOnLineOK( gridLine._length ))
2216       {
2217         ElSLib::Parameters(_cylinder, linCylinder.Point(1) ,_u,_v);
2218         addIntPoint();
2219       }
2220       if ( linCylinder.NbPoints() > 1 )
2221       {
2222         _w = linCylinder.ParamOnConic(2);
2223         if ( isParamOnLineOK( gridLine._length ))
2224         {
2225           ElSLib::Parameters(_cylinder, linCylinder.Point(2) ,_u,_v);
2226           _transition = ( _transition == Trans_OUT ) ? Trans_IN : Trans_OUT;
2227           addIntPoint();
2228         }
2229       }
2230     }
2231   }
2232   //================================================================================
2233   /*
2234    * Intersect a line with a cone
2235    */
2236   void FaceLineIntersector::IntersectWithCone (const GridLine& gridLine)
2237   {
2238     IntAna_IntConicQuad linCone(gridLine._line,_cone);
2239     if ( !linCone.IsDone() ) return;
2240     gp_Pnt P;
2241     gp_Vec du, dv, norm;
2242     for ( int i = 1; i <= linCone.NbPoints(); ++i )
2243     {
2244       _w = linCone.ParamOnConic( i );
2245       if ( !isParamOnLineOK( gridLine._length )) continue;
2246       ElSLib::Parameters(_cone, linCone.Point(i) ,_u,_v);
2247       if ( UVIsOnFace() )
2248       {
2249         ElSLib::D1( _u, _v, _cone, P, du, dv );
2250         norm = du ^ dv;
2251         double normSize2 = norm.SquareMagnitude();
2252         if ( normSize2 > Precision::Angular() * Precision::Angular() )
2253         {
2254           double cos = norm.XYZ() * gridLine._line.Direction().XYZ();
2255           cos /= sqrt( normSize2 );
2256           if ( cos < -Precision::Angular() )
2257             _transition = _transIn;
2258           else if ( cos > Precision::Angular() )
2259             _transition = _transOut;
2260           else
2261             _transition = Trans_TANGENT;
2262         }
2263         else
2264         {
2265           _transition = Trans_APEX;
2266         }
2267         addIntPoint( /*toClassify=*/false);
2268       }
2269     }
2270   }
2271   //================================================================================
2272   /*
2273    * Intersect a line with a sphere
2274    */
2275   void FaceLineIntersector::IntersectWithSphere  (const GridLine& gridLine)
2276   {
2277     IntAna_IntConicQuad linSphere(gridLine._line,_sphere);
2278     if ( linSphere.IsDone() && linSphere.NbPoints() > 0 )
2279     {
2280       _w = linSphere.ParamOnConic(1);
2281       if ( linSphere.NbPoints() == 1 )
2282         _transition = Trans_TANGENT;
2283       else
2284         _transition = _w < linSphere.ParamOnConic(2) ? _transIn : _transOut;
2285       if ( isParamOnLineOK( gridLine._length ))
2286       {
2287         ElSLib::Parameters(_sphere, linSphere.Point(1) ,_u,_v);
2288         addIntPoint();
2289       }
2290       if ( linSphere.NbPoints() > 1 )
2291       {
2292         _w = linSphere.ParamOnConic(2);
2293         if ( isParamOnLineOK( gridLine._length ))
2294         {
2295           ElSLib::Parameters(_sphere, linSphere.Point(2) ,_u,_v);
2296           _transition = ( _transition == Trans_OUT ) ? Trans_IN : Trans_OUT;
2297           addIntPoint();
2298         }
2299       }
2300     }
2301   }
2302   //================================================================================
2303   /*
2304    * Intersect a line with a torus
2305    */
2306   void FaceLineIntersector::IntersectWithTorus   (const GridLine& gridLine)
2307   {
2308     IntAna_IntLinTorus linTorus(gridLine._line,_torus);
2309     if ( !linTorus.IsDone()) return;
2310     gp_Pnt P;
2311     gp_Vec du, dv, norm;
2312     for ( int i = 1; i <= linTorus.NbPoints(); ++i )
2313     {
2314       _w = linTorus.ParamOnLine( i );
2315       if ( !isParamOnLineOK( gridLine._length )) continue;
2316       linTorus.ParamOnTorus( i, _u,_v );
2317       if ( UVIsOnFace() )
2318       {
2319         ElSLib::D1( _u, _v, _torus, P, du, dv );
2320         norm = du ^ dv;
2321         double normSize = norm.Magnitude();
2322         double cos = norm.XYZ() * gridLine._line.Direction().XYZ();
2323         cos /= normSize;
2324         if ( cos < -Precision::Angular() )
2325           _transition = _transIn;
2326         else if ( cos > Precision::Angular() )
2327           _transition = _transOut;
2328         else
2329           _transition = Trans_TANGENT;
2330         addIntPoint( /*toClassify=*/false);
2331       }
2332     }
2333   }
2334   //================================================================================
2335   /*
2336    * Intersect a line with a non-analytical surface
2337    */
2338   void FaceLineIntersector::IntersectWithSurface (const GridLine& gridLine)
2339   {
2340     _surfaceInt->Perform( gridLine._line, 0.0, gridLine._length );
2341     if ( !_surfaceInt->IsDone() ) return;
2342     for ( int i = 1; i <= _surfaceInt->NbPnt(); ++i )
2343     {
2344       _transition = Transition( _surfaceInt->Transition( i ) );
2345       _w = _surfaceInt->WParameter( i );
2346       addIntPoint(/*toClassify=*/false);
2347     }
2348   }
2349   //================================================================================
2350   /*
2351    * check if its face can be safely intersected in a thread
2352    */
2353   bool FaceGridIntersector::IsThreadSafe(set< const Standard_Transient* >& noSafeTShapes) const
2354   {
2355     bool isSafe = true;
2356
2357     // check surface
2358     TopLoc_Location loc;
2359     Handle(Geom_Surface) surf = BRep_Tool::Surface( _face, loc );
2360     Handle(Geom_RectangularTrimmedSurface) ts =
2361       Handle(Geom_RectangularTrimmedSurface)::DownCast( surf );
2362     while( !ts.IsNull() ) {
2363       surf = ts->BasisSurface();
2364       ts = Handle(Geom_RectangularTrimmedSurface)::DownCast(surf);
2365     }
2366     if ( surf->IsKind( STANDARD_TYPE(Geom_BSplineSurface )) ||
2367          surf->IsKind( STANDARD_TYPE(Geom_BezierSurface )))
2368       if ( !noSafeTShapes.insert( _face.TShape().get() ).second )
2369         isSafe = false;
2370
2371     double f, l;
2372     TopExp_Explorer exp( _face, TopAbs_EDGE );
2373     for ( ; exp.More(); exp.Next() )
2374     {
2375       bool edgeIsSafe = true;
2376       const TopoDS_Edge& e = TopoDS::Edge( exp.Current() );
2377       // check 3d curve
2378       {
2379         Handle(Geom_Curve) c = BRep_Tool::Curve( e, loc, f, l);
2380         if ( !c.IsNull() )
2381         {
2382           Handle(Geom_TrimmedCurve) tc = Handle(Geom_TrimmedCurve)::DownCast(c);
2383           while( !tc.IsNull() ) {
2384             c = tc->BasisCurve();
2385             tc = Handle(Geom_TrimmedCurve)::DownCast(c);
2386           }
2387           if ( c->IsKind( STANDARD_TYPE(Geom_BSplineCurve )) ||
2388                c->IsKind( STANDARD_TYPE(Geom_BezierCurve )))
2389             edgeIsSafe = false;
2390         }
2391       }
2392       // check 2d curve
2393       if ( edgeIsSafe )
2394       {
2395         Handle(Geom2d_Curve) c2 = BRep_Tool::CurveOnSurface( e, surf, loc, f, l);
2396         if ( !c2.IsNull() )
2397         {
2398           Handle(Geom2d_TrimmedCurve) tc = Handle(Geom2d_TrimmedCurve)::DownCast(c2);
2399           while( !tc.IsNull() ) {
2400             c2 = tc->BasisCurve();
2401             tc = Handle(Geom2d_TrimmedCurve)::DownCast(c2);
2402           }
2403           if ( c2->IsKind( STANDARD_TYPE(Geom2d_BSplineCurve )) ||
2404                c2->IsKind( STANDARD_TYPE(Geom2d_BezierCurve )))
2405             edgeIsSafe = false;
2406         }
2407       }
2408       if ( !edgeIsSafe && !noSafeTShapes.insert( e.TShape().get() ).second )
2409         isSafe = false;
2410     }
2411     return isSafe;
2412   }
2413   //================================================================================
2414   /*!
2415    * \brief Creates topology of the hexahedron
2416    */
2417   Hexahedron::Hexahedron(Grid* grid)
2418     : _grid( grid ), _nbFaceIntNodes(0), _hasTooSmall( false )
2419   {
2420     _polygons.reserve(100); // to avoid reallocation;
2421
2422     //set nodes shift within grid->_nodes from the node 000 
2423     size_t dx = _grid->NodeIndexDX();
2424     size_t dy = _grid->NodeIndexDY();
2425     size_t dz = _grid->NodeIndexDZ();
2426     size_t i000 = 0;
2427     size_t i100 = i000 + dx;
2428     size_t i010 = i000 + dy;
2429     size_t i110 = i010 + dx;
2430     size_t i001 = i000 + dz;
2431     size_t i101 = i100 + dz;
2432     size_t i011 = i010 + dz;
2433     size_t i111 = i110 + dz;
2434     grid->_nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V000 )] = i000;
2435     grid->_nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V100 )] = i100;
2436     grid->_nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V010 )] = i010;
2437     grid->_nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V110 )] = i110;
2438     grid->_nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V001 )] = i001;
2439     grid->_nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V101 )] = i101;
2440     grid->_nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V011 )] = i011;
2441     grid->_nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V111 )] = i111;
2442
2443     vector< int > idVec;
2444     // set nodes to links
2445     for ( int linkID = SMESH_Block::ID_Ex00; linkID <= SMESH_Block::ID_E11z; ++linkID )
2446     {
2447       SMESH_Block::GetEdgeVertexIDs( linkID, idVec );
2448       _Link& link = _hexLinks[ SMESH_Block::ShapeIndex( linkID )];
2449       link._nodes[0] = &_hexNodes[ SMESH_Block::ShapeIndex( idVec[0] )];
2450       link._nodes[1] = &_hexNodes[ SMESH_Block::ShapeIndex( idVec[1] )];
2451     }
2452
2453     // set links to faces
2454     int interlace[4] = { 0, 3, 1, 2 }; // to walk by links around a face: { u0, 1v, u1, 0v }
2455     for ( int faceID = SMESH_Block::ID_Fxy0; faceID <= SMESH_Block::ID_F1yz; ++faceID )
2456     {
2457       _Face& quad = _hexQuads[ SMESH_Block::ShapeIndex( faceID )];
2458       quad._name = (SMESH_Block::TShapeID) faceID;
2459
2460       SMESH_Block::GetFaceEdgesIDs( faceID, idVec );
2461       bool revFace = ( faceID == SMESH_Block::ID_Fxy0 ||
2462                        faceID == SMESH_Block::ID_Fx1z ||
2463                        faceID == SMESH_Block::ID_F0yz );
2464       quad._links.resize(4);
2465       vector<_OrientedLink>::iterator         frwLinkIt = quad._links.begin();
2466       vector<_OrientedLink>::reverse_iterator revLinkIt = quad._links.rbegin();
2467       for ( int i = 0; i < 4; ++i )
2468       {
2469         bool revLink = revFace;
2470         if ( i > 1 ) // reverse links u1 and v0
2471           revLink = !revLink;
2472         _OrientedLink& link = revFace ? *revLinkIt++ : *frwLinkIt++;
2473         link = _OrientedLink( & _hexLinks[ SMESH_Block::ShapeIndex( idVec[interlace[i]] )],
2474                               revLink );
2475       }
2476     }
2477   }
2478   //================================================================================
2479   /*!
2480    * \brief Copy constructor
2481    */
2482   Hexahedron::Hexahedron( const Hexahedron& other, size_t i, size_t j, size_t k, int cellID )
2483     :_grid( other._grid ), _nbFaceIntNodes(0), _i( i ), _j( j ), _k( k ), _hasTooSmall( false )
2484   {
2485     _polygons.reserve(100); // to avoid reallocation;
2486
2487     // copy topology
2488     for ( int i = 0; i < 12; ++i )
2489     {
2490       const _Link& srcLink = other._hexLinks[ i ];
2491       _Link&       tgtLink = this->_hexLinks[ i ];
2492       tgtLink._nodes[0] = _hexNodes + ( srcLink._nodes[0] - other._hexNodes );
2493       tgtLink._nodes[1] = _hexNodes + ( srcLink._nodes[1] - other._hexNodes );
2494     }
2495
2496     for ( int i = 0; i < 6; ++i )
2497     {
2498       const _Face& srcQuad = other._hexQuads[ i ];
2499       _Face&       tgtQuad = this->_hexQuads[ i ];
2500       tgtQuad._name = srcQuad._name;
2501       tgtQuad._links.resize(4);
2502       for ( int j = 0; j < 4; ++j )
2503       {
2504         const _OrientedLink& srcLink = srcQuad._links[ j ];
2505         _OrientedLink&       tgtLink = tgtQuad._links[ j ];
2506         tgtLink._reverse = srcLink._reverse;
2507         tgtLink._link    = _hexLinks + ( srcLink._link - other._hexLinks );
2508       }
2509     }
2510     
2511     if (SALOME::VerbosityActivated())
2512       _cellID = cellID;
2513   }
2514
2515   //================================================================================
2516   /*!
2517    * \brief Return IDs of SOLIDs interfering with this Hexahedron
2518    */
2519   size_t Hexahedron::getSolids( TGeomID ids[] )
2520   {
2521     if ( _grid->_geometry.IsOneSolid() )
2522     {
2523       ids[0] = _grid->GetSolid()->ID();
2524       return 1;
2525     }
2526     // count intersection points belonging to each SOLID
2527     TID2Nb id2NbPoints;
2528     id2NbPoints.reserve( 3 );
2529
2530     _origNodeInd = _grid->NodeIndex( _i,_j,_k );
2531     for ( int iN = 0; iN < 8; ++iN )
2532     {
2533       _hexNodes[iN]._node     = _grid->_nodes   [ _origNodeInd + _grid->_nodeShift[iN] ];
2534       _hexNodes[iN]._intPoint = _grid->_gridIntP[ _origNodeInd + _grid->_nodeShift[iN] ];
2535
2536       if ( _hexNodes[iN]._intPoint ) // intersection with a FACE
2537       {
2538         for ( size_t iF = 0; iF < _hexNodes[iN]._intPoint->_faceIDs.size(); ++iF )
2539         {
2540           const vector< TGeomID > & solidIDs =
2541             _grid->GetSolidIDs( _hexNodes[iN]._intPoint->_faceIDs[iF] );
2542           for ( size_t i = 0; i < solidIDs.size(); ++i )
2543             insertAndIncrement( solidIDs[i], id2NbPoints );
2544         }
2545       }
2546       else if ( _hexNodes[iN]._node ) // node inside a SOLID
2547       {
2548         insertAndIncrement( _hexNodes[iN]._node->GetShapeID(), id2NbPoints );
2549       }
2550     }
2551
2552     for ( int iL = 0; iL < 12; ++iL )
2553     {
2554       const _Link& link = _hexLinks[ iL ];
2555       for ( size_t iP = 0; iP < link._fIntPoints.size(); ++iP )
2556       {
2557         for ( size_t iF = 0; iF < link._fIntPoints[iP]->_faceIDs.size(); ++iF )
2558         {
2559           const vector< TGeomID > & solidIDs =
2560             _grid->GetSolidIDs( link._fIntPoints[iP]->_faceIDs[iF] );
2561           for ( size_t i = 0; i < solidIDs.size(); ++i )
2562             insertAndIncrement( solidIDs[i], id2NbPoints );
2563         }
2564       }
2565     }
2566
2567     for ( size_t iP = 0; iP < _eIntPoints.size(); ++iP )
2568     {
2569       const vector< TGeomID > & solidIDs = _grid->GetSolidIDs( _eIntPoints[iP]->_shapeID );
2570       for ( size_t i = 0; i < solidIDs.size(); ++i )
2571         insertAndIncrement( solidIDs[i], id2NbPoints );
2572     }
2573
2574     size_t nbSolids = 0;
2575     for ( TID2Nb::iterator id2nb = id2NbPoints.begin(); id2nb != id2NbPoints.end(); ++id2nb )
2576       if ( id2nb->second >= 3 )
2577         ids[ nbSolids++ ] = id2nb->first;
2578
2579     return nbSolids;
2580   }
2581
2582   //================================================================================
2583   /*!
2584    * \brief Count cuts by INTERNAL FACEs and set _Node::_isInternalFlags
2585    */
2586   bool Hexahedron::isCutByInternalFace( IsInternalFlag & maxFlag )
2587   {
2588     TID2Nb id2NbPoints;
2589     id2NbPoints.reserve( 3 );
2590
2591     for ( size_t iN = 0; iN < _intNodes.size(); ++iN )
2592       for ( size_t iF = 0; iF < _intNodes[iN]._intPoint->_faceIDs.size(); ++iF )
2593       {
2594         if ( _grid->IsInternal( _intNodes[iN]._intPoint->_faceIDs[iF]))
2595           insertAndIncrement( _intNodes[iN]._intPoint->_faceIDs[iF], id2NbPoints );
2596       }
2597     for ( size_t iN = 0; iN < 8; ++iN )
2598       if ( _hexNodes[iN]._intPoint )
2599         for ( size_t iF = 0; iF < _hexNodes[iN]._intPoint->_faceIDs.size(); ++iF )
2600         {
2601           if ( _grid->IsInternal( _hexNodes[iN]._intPoint->_faceIDs[iF]))
2602             insertAndIncrement( _hexNodes[iN]._intPoint->_faceIDs[iF], id2NbPoints );
2603         }
2604
2605     maxFlag = IS_NOT_INTERNAL;
2606     for ( TID2Nb::iterator id2nb = id2NbPoints.begin(); id2nb != id2NbPoints.end(); ++id2nb )
2607     {
2608       TGeomID        intFace = id2nb->first;
2609       IsInternalFlag intFlag = ( id2nb->second >= 3 ? IS_CUT_BY_INTERNAL_FACE : IS_INTERNAL );
2610       if ( intFlag > maxFlag )
2611         maxFlag = intFlag;
2612
2613       for ( size_t iN = 0; iN < _intNodes.size(); ++iN )
2614         if ( _intNodes[iN].IsOnFace( intFace ))
2615           _intNodes[iN].SetInternal( intFlag );
2616
2617       for ( size_t iN = 0; iN < 8; ++iN )
2618         if ( _hexNodes[iN].IsOnFace( intFace ))
2619           _hexNodes[iN].SetInternal( intFlag );
2620     }
2621
2622     return maxFlag;
2623   }
2624
2625   //================================================================================
2626   /*!
2627    * \brief Return any FACE interfering with this Hexahedron
2628    */
2629   TGeomID Hexahedron::getAnyFace() const
2630   {
2631     TID2Nb id2NbPoints;
2632     id2NbPoints.reserve( 3 );
2633
2634     for ( size_t iN = 0; iN < _intNodes.size(); ++iN )
2635       for ( size_t iF = 0; iF < _intNodes[iN]._intPoint->_faceIDs.size(); ++iF )
2636         insertAndIncrement( _intNodes[iN]._intPoint->_faceIDs[iF], id2NbPoints );
2637
2638     for ( size_t iN = 0; iN < 8; ++iN )
2639       if ( _hexNodes[iN]._intPoint )
2640         for ( size_t iF = 0; iF < _hexNodes[iN]._intPoint->_faceIDs.size(); ++iF )
2641           insertAndIncrement( _hexNodes[iN]._intPoint->_faceIDs[iF], id2NbPoints );
2642
2643     for ( unsigned int minNb = 3; minNb > 0; --minNb )
2644       for ( TID2Nb::iterator id2nb = id2NbPoints.begin(); id2nb != id2NbPoints.end(); ++id2nb )
2645         if ( id2nb->second >= minNb )
2646           return id2nb->first;
2647
2648     return 0;
2649   }
2650
2651   //================================================================================
2652   /*!
2653    * \brief Initializes IJK by Hexahedron index
2654    */
2655   void Hexahedron::setIJK( size_t iCell )
2656   {
2657     size_t iNbCell = _grid->_coords[0].size() - 1;
2658     size_t jNbCell = _grid->_coords[1].size() - 1;
2659     _i = iCell % iNbCell;
2660     _j = ( iCell % ( iNbCell * jNbCell )) / iNbCell;
2661     _k = iCell / iNbCell / jNbCell;
2662   }
2663
2664   //================================================================================
2665   /*!
2666    * \brief Initializes its data by given grid cell (countered from zero)
2667    */
2668   void Hexahedron::init( size_t iCell )
2669   {
2670     setIJK( iCell );
2671     init( _i, _j, _k );
2672   }
2673
2674   //================================================================================
2675   /*!
2676    * \brief Initializes its data by given grid cell nodes and intersections
2677    */
2678   void Hexahedron::init( size_t i, size_t j, size_t k, const Solid* solid )
2679   {
2680     _i = i; _j = j; _k = k;
2681
2682     bool isCompute = solid;
2683     if ( !solid )
2684       solid = _grid->GetSolid();
2685
2686     // set nodes of grid to nodes of the hexahedron and
2687     // count nodes at hexahedron corners located IN and ON geometry
2688     _nbCornerNodes = _nbBndNodes = 0;
2689     _origNodeInd   = _grid->NodeIndex( i,j,k );
2690     for ( int iN = 0; iN < 8; ++iN )
2691     {
2692       _hexNodes[iN]._isInternalFlags = 0;
2693
2694       // Grid  node 
2695       _hexNodes[iN]._node     = _grid->_nodes   [ _origNodeInd + _grid->_nodeShift[iN] ];
2696       _hexNodes[iN]._intPoint = _grid->_gridIntP[ _origNodeInd + _grid->_nodeShift[iN] ];
2697
2698       if ( _grid->_allBorderNodes[ _origNodeInd + _grid->_nodeShift[iN] ] ) 
2699         _hexNodes[iN]._boundaryCornerNode = _grid->_allBorderNodes [ _origNodeInd + _grid->_nodeShift[iN] ];
2700       
2701       if ( _hexNodes[iN]._node && !solid->Contains( _hexNodes[iN]._node->GetShapeID() ))
2702         _hexNodes[iN]._node = 0;
2703
2704       if ( _hexNodes[iN]._intPoint && !solid->ContainsAny( _hexNodes[iN]._intPoint->_faceIDs ))
2705         _hexNodes[iN]._intPoint = 0;
2706
2707       _nbCornerNodes += bool( _hexNodes[iN]._node );
2708       _nbBndNodes    += bool( _hexNodes[iN]._intPoint );
2709     }
2710     _sideLength[0] = _grid->_coords[0][i+1] - _grid->_coords[0][i];
2711     _sideLength[1] = _grid->_coords[1][j+1] - _grid->_coords[1][j];
2712     _sideLength[2] = _grid->_coords[2][k+1] - _grid->_coords[2][k];
2713
2714     _intNodes.clear();
2715     _vIntNodes.clear();
2716
2717     if ( !isCompute )
2718       return;
2719
2720     if ( _nbFaceIntNodes + _eIntPoints.size()                  > 0 &&
2721          _nbFaceIntNodes + _eIntPoints.size() + _nbCornerNodes > 3)
2722     {
2723       _intNodes.reserve( 3 * ( _nbBndNodes + _nbFaceIntNodes + _eIntPoints.size() ));
2724
2725       // this method can be called in parallel, so use own helper
2726       SMESH_MesherHelper helper( *_grid->_helper->GetMesh() );
2727
2728       // Create sub-links (_Link::_splits) by splitting links with _Link::_fIntPoints
2729       // ---------------------------------------------------------------
2730       _Link split;
2731       for ( int iLink = 0; iLink < 12; ++iLink )
2732       {
2733         _Link& link = _hexLinks[ iLink ];
2734         link._fIntNodes.clear();
2735         link._fIntNodes.reserve( link._fIntPoints.size() );
2736         for ( size_t i = 0; i < link._fIntPoints.size(); ++i )
2737           if ( solid->ContainsAny( link._fIntPoints[i]->_faceIDs ))
2738           {
2739             _intNodes.push_back( _Node( 0, link._fIntPoints[i] ));
2740             link._fIntNodes.push_back( & _intNodes.back() );
2741           }
2742
2743         link._splits.clear();
2744         split._nodes[ 0 ] = link._nodes[0];
2745         bool isOut = ( ! link._nodes[0]->Node() );
2746         bool checkTransition;
2747         for ( size_t i = 0; i < link._fIntNodes.size(); ++i )
2748         {
2749           const bool isGridNode = ( ! link._fIntNodes[i]->Node() );
2750           if ( !isGridNode ) // intersection non-coincident with a grid node
2751           {
2752             if ( split._nodes[ 0 ]->Node() && !isOut )
2753             {
2754               split._nodes[ 1 ] = link._fIntNodes[i];
2755               link._splits.push_back( split );
2756             }
2757             split._nodes[ 0 ] = link._fIntNodes[i];
2758             checkTransition = true;
2759           }
2760           else // FACE intersection coincident with a grid node (at link ends)
2761           {
2762             checkTransition = ( i == 0 && link._nodes[0]->Node() );
2763           }
2764           if ( checkTransition )
2765           {
2766             const vector< TGeomID >& faceIDs = link._fIntNodes[i]->_intPoint->_faceIDs;
2767             if ( _grid->IsInternal( faceIDs.back() ))
2768               isOut = false;
2769             else if ( faceIDs.size() > 1 || _eIntPoints.size() > 0 )
2770               isOut = isOutPoint( link, i, helper, solid );
2771             else
2772             {
2773               bool okTransi = _grid->IsCorrectTransition( faceIDs[0], solid );
2774               switch ( link._fIntNodes[i]->FaceIntPnt()->_transition ) {
2775               case Trans_OUT: isOut = okTransi;  break;
2776               case Trans_IN : isOut = !okTransi; break;
2777               default:
2778                 isOut = isOutPoint( link, i, helper, solid );
2779               }
2780             }
2781           }
2782         }
2783         if ( link._nodes[ 1 ]->Node() && split._nodes[ 0 ]->Node() && !isOut )
2784         {
2785           split._nodes[ 1 ] = link._nodes[1];
2786           link._splits.push_back( split );
2787         }
2788       }
2789
2790       // Create _Node's at intersections with EDGEs.
2791       // --------------------------------------------
2792       // 1) add this->_eIntPoints to _Face::_eIntNodes
2793       // 2) fill _intNodes and _vIntNodes
2794       //
2795       const double tol2 = _grid->_tol * _grid->_tol * 4;
2796       int facets[3], nbFacets, subEntity;
2797
2798       for ( int iF = 0; iF < 6; ++iF )
2799         _hexQuads[ iF ]._eIntNodes.clear();
2800
2801       for ( size_t iP = 0; iP < _eIntPoints.size(); ++iP )
2802       {
2803         if ( !solid->ContainsAny( _eIntPoints[iP]->_faceIDs ))
2804           continue;
2805         nbFacets = getEntity( _eIntPoints[iP], facets, subEntity );
2806         _Node* equalNode = 0;
2807         switch( nbFacets ) {
2808         case 1: // in a _Face
2809         {
2810           _Face& quad = _hexQuads[ facets[0] - SMESH_Block::ID_FirstF ];
2811           equalNode = findEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
2812           if ( equalNode ) {
2813             equalNode->Add( _eIntPoints[ iP ] );
2814           }
2815           else {
2816             _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
2817             quad._eIntNodes.push_back( & _intNodes.back() );
2818           }
2819           break;
2820         }
2821         case 2: // on a _Link
2822         {
2823           _Link& link = _hexLinks[ subEntity - SMESH_Block::ID_FirstE ];
2824           if ( link._splits.size() > 0 )
2825           {
2826             equalNode = findEqualNode( link._fIntNodes, _eIntPoints[ iP ], tol2 );
2827             if ( equalNode )
2828               equalNode->Add( _eIntPoints[ iP ] );
2829             else if ( link._splits.size() == 1 &&
2830                       link._splits[0]._nodes[0] &&
2831                       link._splits[0]._nodes[1] )
2832               link._splits.clear(); // hex edge is divided by _eIntPoints[iP]
2833           }
2834           //else
2835           if ( !equalNode )
2836           {
2837             _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
2838             bool newNodeUsed = false;
2839             for ( int iF = 0; iF < 2; ++iF )
2840             {
2841               _Face& quad = _hexQuads[ facets[iF] - SMESH_Block::ID_FirstF ];
2842               equalNode = findEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
2843               if ( equalNode ) {
2844                 equalNode->Add( _eIntPoints[ iP ] );
2845               }
2846               else {
2847                 quad._eIntNodes.push_back( & _intNodes.back() );
2848                 newNodeUsed = true;
2849               }
2850             }
2851             if ( !newNodeUsed )
2852               _intNodes.pop_back();
2853           }
2854           break;
2855         }
2856         case 3: // at a corner
2857         {
2858           _Node& node = _hexNodes[ subEntity - SMESH_Block::ID_FirstV ];
2859           if ( node.Node() )
2860           {
2861             if ( node._intPoint )
2862               node._intPoint->Add( _eIntPoints[ iP ]->_faceIDs, _eIntPoints[ iP ]->_node );
2863           }
2864           else
2865           {
2866             _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
2867             for ( int iF = 0; iF < 3; ++iF )
2868             {
2869               _Face& quad = _hexQuads[ facets[iF] - SMESH_Block::ID_FirstF ];
2870               equalNode = findEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
2871               if ( equalNode ) {
2872                 equalNode->Add( _eIntPoints[ iP ] );
2873               }
2874               else {
2875                 quad._eIntNodes.push_back( & _intNodes.back() );
2876               }
2877             }
2878           }
2879           break;
2880         }
2881         } // switch( nbFacets )
2882
2883         if ( nbFacets == 0 ||
2884              _grid->ShapeType( _eIntPoints[ iP ]->_shapeID ) == TopAbs_VERTEX )
2885         {
2886           equalNode = findEqualNode( _vIntNodes, _eIntPoints[ iP ], tol2 );
2887           if ( equalNode ) {
2888             equalNode->Add( _eIntPoints[ iP ] );
2889           }
2890           else if ( nbFacets == 0 ) {
2891             if ( _intNodes.empty() || _intNodes.back().EdgeIntPnt() != _eIntPoints[ iP ])
2892               _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
2893             _vIntNodes.push_back( & _intNodes.back() );
2894           }
2895         }
2896       } // loop on _eIntPoints
2897     }
2898
2899     else if (( 3 < _nbCornerNodes && _nbCornerNodes < 8 ) || // _nbFaceIntNodes == 0
2900              ( !_grid->_geometry.IsOneSolid() ))
2901     {
2902       _Link split;
2903       // create sub-links (_splits) of whole links
2904       for ( int iLink = 0; iLink < 12; ++iLink )
2905       {
2906         _Link& link = _hexLinks[ iLink ];
2907         link._splits.clear();
2908         if ( link._nodes[ 0 ]->Node() && link._nodes[ 1 ]->Node() )
2909         {
2910           split._nodes[ 0 ] = link._nodes[0];
2911           split._nodes[ 1 ] = link._nodes[1];
2912           link._splits.push_back( split );
2913         }
2914       }
2915     }
2916     return;
2917
2918   } // init( _i, _j, _k )
2919
2920   //================================================================================
2921   /*!
2922    * \brief Compute mesh volumes resulted from intersection of the Hexahedron
2923    */
2924   void Hexahedron::computeElements( const Solid* solid, int solidIndex )
2925   {
2926     if ( !solid )
2927     {
2928       solid = _grid->GetSolid();
2929       if ( !_grid->_geometry.IsOneSolid() )
2930       {
2931         TGeomID solidIDs[20] = { 0 };
2932         size_t nbSolids = getSolids( solidIDs );
2933         if ( nbSolids > 1 )
2934         {
2935           for ( size_t i = 0; i < nbSolids; ++i )
2936           {
2937             solid = _grid->GetSolid( solidIDs[i] );
2938             computeElements( solid, i );
2939             if ( !_volumeDefs._nodes.empty() && i < nbSolids - 1 )
2940               _volumeDefs.SetNext( new _volumeDef( _volumeDefs ));
2941           }
2942           return;
2943         }
2944         solid = _grid->GetSolid( solidIDs[0] );
2945       }
2946     }
2947
2948     init( _i, _j, _k, solid ); // get nodes and intersections from grid nodes and split links
2949
2950     int nbIntersections = _nbFaceIntNodes + _eIntPoints.size();
2951     if ( _nbCornerNodes + nbIntersections < 4 )
2952       return;
2953
2954     if ( _nbBndNodes == _nbCornerNodes && nbIntersections == 0 && isInHole() )
2955       return; // cell is in a hole
2956
2957     IsInternalFlag intFlag = IS_NOT_INTERNAL;
2958     if ( solid->HasInternalFaces() && this->isCutByInternalFace( intFlag ))
2959     {
2960       for ( _SplitIterator it( _hexLinks ); it.More(); it.Next() )
2961       {
2962         if ( compute( solid, intFlag ))
2963           _volumeDefs.SetNext( new _volumeDef( _volumeDefs ));
2964       }
2965     }
2966     else
2967     {
2968       if ( solidIndex >= 0 )
2969         intFlag = IS_CUT_BY_INTERNAL_FACE;
2970
2971       compute( solid, intFlag );
2972     }
2973   }
2974
2975   //================================================================================
2976   /*!
2977    * \brief Compute mesh volumes resulted from intersection of the Hexahedron
2978    */
2979   bool Hexahedron::compute( const Solid* solid, const IsInternalFlag intFlag )
2980   {
2981     _polygons.clear();
2982     _polygons.reserve( 20 );
2983
2984     for ( int iN = 0; iN < 8; ++iN )
2985       _hexNodes[iN]._usedInFace = 0;
2986
2987     if ( intFlag & IS_CUT_BY_INTERNAL_FACE && !_grid->_toAddEdges ) // Issue #19913
2988       preventVolumesOverlapping();
2989
2990     std::set< TGeomID > concaveFaces; // to avoid connecting nodes laying on them
2991
2992     if ( solid->HasConcaveVertex() )
2993     {
2994       for ( const E_IntersectPoint* ip : _eIntPoints )
2995       {
2996         if ( const ConcaveFace* cf = solid->GetConcave( ip->_shapeID ))
2997           if ( this->hasEdgesAround( cf ))
2998             concaveFaces.insert( cf->_concaveFace );
2999       }
3000       if ( concaveFaces.empty() || concaveFaces.size() * 3  < _eIntPoints.size() )
3001         for ( const _Node& hexNode: _hexNodes )
3002         {
3003           if ( hexNode._node && hexNode._intPoint && hexNode._intPoint->_faceIDs.size() >= 3 )
3004             if ( const ConcaveFace* cf = solid->GetConcave( hexNode._node->GetShapeID() ))
3005               if ( this->hasEdgesAround( cf ))
3006                 concaveFaces.insert( cf->_concaveFace );
3007         }
3008     }
3009
3010     // Create polygons from quadrangles
3011     // --------------------------------
3012
3013     vector< _OrientedLink > splits;
3014     vector<_Node*>          chainNodes;
3015     _Face*                  coplanarPolyg;
3016
3017     const bool hasEdgeIntersections = !_eIntPoints.empty();
3018     const bool toCheckSideDivision = isImplementEdges() || intFlag & IS_CUT_BY_INTERNAL_FACE;
3019
3020     for ( int iF = 0; iF < 6; ++iF ) // loop on 6 sides of a hexahedron
3021     {
3022       _Face& quad = _hexQuads[ iF ] ;
3023
3024       _polygons.resize( _polygons.size() + 1 );
3025       _Face* polygon = &_polygons.back();
3026       polygon->_polyLinks.reserve( 20 );
3027       polygon->_name = quad._name;
3028
3029       splits.clear();
3030       for ( int iE = 0; iE < 4; ++iE ) // loop on 4 sides of a quadrangle
3031         for ( size_t iS = 0; iS < quad._links[ iE ].NbResultLinks(); ++iS )
3032           splits.push_back( quad._links[ iE ].ResultLink( iS ));
3033
3034       if ( splits.size() == 4 &&
3035            isQuadOnFace( iF )) // check if a quad on FACE is not split
3036       {
3037         polygon->_links.swap( splits );
3038         continue; // goto the next quad
3039       }
3040
3041       // add splits of links to a polygon and add _polyLinks to make
3042       // polygon's boundary closed
3043
3044       int nbSplits = splits.size();
3045       if (( nbSplits == 1 ) &&
3046           ( quad._eIntNodes.empty() ||
3047             splits[0].FirstNode()->IsLinked( splits[0].LastNode()->_intPoint )))
3048         //( quad._eIntNodes.empty() || _nbCornerNodes + nbIntersections > 6 ))
3049         nbSplits = 0;
3050
3051       for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
3052         if ( quad._eIntNodes[ iP ]->IsUsedInFace( polygon ))
3053           quad._eIntNodes[ iP ]->_usedInFace = 0;
3054
3055       size_t nbUsedEdgeNodes = 0;
3056       _Face* prevPolyg = 0; // polygon previously created from this quad
3057
3058       while ( nbSplits > 0 )
3059       {
3060         size_t iS = 0;
3061         while ( !splits[ iS ] )
3062           ++iS;
3063
3064         if ( !polygon->_links.empty() )
3065         {
3066           _polygons.resize( _polygons.size() + 1 );
3067           polygon = &_polygons.back();
3068           polygon->_polyLinks.reserve( 20 );
3069           polygon->_name = quad._name;
3070         }
3071         polygon->_links.push_back( splits[ iS ] );
3072         splits[ iS++ ]._link = 0;
3073         --nbSplits;
3074
3075         _Node* nFirst = polygon->_links.back().FirstNode();
3076         _Node *n1,*n2 = polygon->_links.back().LastNode();
3077         for ( ; nFirst != n2 && iS < splits.size(); ++iS )
3078         {
3079           _OrientedLink& split = splits[ iS ];
3080           if ( !split ) continue;
3081
3082           n1 = split.FirstNode();
3083           if ( n1 == n2 &&
3084                n1->_intPoint &&
3085                (( n1->_intPoint->_faceIDs.size() > 1 && toCheckSideDivision ) ||
3086                 ( n1->_isInternalFlags )))
3087           {
3088             // n1 is at intersection with EDGE
3089             if ( findChainOnEdge( splits, polygon->_links.back(), split, concaveFaces,
3090                                   iS, quad, chainNodes ))
3091             {
3092               for ( size_t i = 1; i < chainNodes.size(); ++i )
3093                 polygon->AddPolyLink( chainNodes[i-1], chainNodes[i], prevPolyg );
3094               if ( chainNodes.back() != n1 ) // not a partial cut by INTERNAL FACE
3095               {
3096                 prevPolyg = polygon;
3097                 n2 = chainNodes.back();
3098                 continue;
3099               }
3100             }
3101           }
3102           else if ( n1 != n2 )
3103           {
3104             // try to connect to intersections with EDGEs
3105             if ( quad._eIntNodes.size() > nbUsedEdgeNodes  &&
3106                  findChain( n2, n1, quad, chainNodes ))
3107             {
3108               for ( size_t i = 1; i < chainNodes.size(); ++i )
3109               {
3110                 polygon->AddPolyLink( chainNodes[i-1], chainNodes[i] );
3111                 nbUsedEdgeNodes += ( chainNodes[i]->IsUsedInFace( polygon ));
3112               }
3113               if ( chainNodes.back() != n1 )
3114               {
3115                 n2 = chainNodes.back();
3116                 --iS;
3117                 continue;
3118               }
3119             }
3120             // try to connect to a split ending on the same FACE
3121             else
3122             {
3123               _OrientedLink foundSplit;
3124               for ( size_t i = iS; i < splits.size() && !foundSplit; ++i )
3125                 if (( foundSplit = splits[ i ]) &&
3126                     ( n2->IsLinked( foundSplit.FirstNode()->_intPoint )))
3127                 {
3128                   iS = i - 1;
3129                 }
3130                 else
3131                 {
3132                   foundSplit._link = 0;
3133                 }
3134               if ( foundSplit )
3135               {
3136                 if ( n2 != foundSplit.FirstNode() )
3137                 {
3138                   polygon->AddPolyLink( n2, foundSplit.FirstNode() );
3139                   n2 = foundSplit.FirstNode();
3140                 }
3141                 continue;
3142               }
3143               else
3144               {
3145                 if ( n2->IsLinked( nFirst->_intPoint ))
3146                   break;
3147                 polygon->AddPolyLink( n2, n1, prevPolyg );
3148               }
3149             }
3150           } // if ( n1 != n2 )
3151
3152           polygon->_links.push_back( split );
3153           split._link = 0;
3154           --nbSplits;
3155           n2 = polygon->_links.back().LastNode();
3156
3157         } // loop on splits
3158
3159         if ( nFirst != n2 ) // close a polygon
3160         {
3161           if ( !findChain( n2, nFirst, quad, chainNodes ))
3162           {
3163             if ( !closePolygon( polygon, chainNodes ))
3164               if ( !isImplementEdges() )
3165                 chainNodes.push_back( nFirst );
3166           }
3167           for ( size_t i = 1; i < chainNodes.size(); ++i )
3168           {
3169             polygon->AddPolyLink( chainNodes[i-1], chainNodes[i], prevPolyg );
3170             nbUsedEdgeNodes += bool( chainNodes[i]->IsUsedInFace( polygon ));
3171           }
3172         }
3173
3174         if ( polygon->_links.size() < 3 && nbSplits > 0 )
3175         {
3176           polygon->_polyLinks.clear();
3177           polygon->_links.clear();
3178         }
3179       } // while ( nbSplits > 0 )
3180
3181       if ( polygon->_links.size() < 3 )
3182       {
3183         _polygons.pop_back();
3184       }
3185     }  // loop on 6 hexahedron sides
3186
3187     // Create polygons closing holes in a polyhedron
3188     // ----------------------------------------------
3189
3190     // clear _usedInFace
3191     for ( size_t iN = 0; iN < _intNodes.size(); ++iN )
3192       _intNodes[ iN ]._usedInFace = 0;
3193
3194     // add polygons to their links and mark used nodes
3195     for ( size_t iP = 0; iP < _polygons.size(); ++iP )
3196     {
3197       _Face& polygon = _polygons[ iP ];
3198       for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
3199       {
3200         polygon._links[ iL ].AddFace( &polygon );
3201         polygon._links[ iL ].FirstNode()->_usedInFace = &polygon;
3202       }
3203     }
3204     // find free links
3205     vector< _OrientedLink* > freeLinks;
3206     freeLinks.reserve(20);
3207     for ( size_t iP = 0; iP < _polygons.size(); ++iP )
3208     {
3209       _Face& polygon = _polygons[ iP ];
3210       for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
3211         if ( polygon._links[ iL ].NbFaces() < 2 )
3212           freeLinks.push_back( & polygon._links[ iL ]);
3213     }
3214     int nbFreeLinks = freeLinks.size();
3215     if ( nbFreeLinks == 1 ) return false;
3216
3217     // put not used intersection nodes to _vIntNodes
3218     int nbVertexNodes = 0; // nb not used vertex nodes
3219     {
3220       for ( size_t iN = 0; iN < _vIntNodes.size(); ++iN )
3221         nbVertexNodes += ( !_vIntNodes[ iN ]->IsUsedInFace() );
3222
3223       const double tol = 1e-3 * Min( Min( _sideLength[0], _sideLength[1] ), _sideLength[0] );
3224       for ( size_t iN = _nbFaceIntNodes; iN < _intNodes.size(); ++iN )
3225       {
3226         if ( _intNodes[ iN ].IsUsedInFace() ) continue;
3227         if ( dynamic_cast< const F_IntersectPoint* >( _intNodes[ iN ]._intPoint )) continue;
3228         _Node* equalNode =
3229           findEqualNode( _vIntNodes, _intNodes[ iN ].EdgeIntPnt(), tol*tol );
3230         if ( !equalNode )
3231         {
3232           _vIntNodes.push_back( &_intNodes[ iN ]);
3233           ++nbVertexNodes;
3234         }
3235       }
3236     }
3237
3238     std::set<TGeomID> usedFaceIDs;
3239     std::vector< TGeomID > faces;
3240     TGeomID curFace = 0;
3241     const size_t nbQuadPolygons = _polygons.size();
3242     E_IntersectPoint ipTmp;
3243     std::map< TGeomID, std::vector< const B_IntersectPoint* > > tmpAddedFace; // face added to _intPoint
3244
3245     // create polygons by making closed chains of free links
3246     size_t iPolygon = _polygons.size();
3247     while ( nbFreeLinks > 0 )
3248     {
3249       if ( iPolygon == _polygons.size() )
3250       {
3251         _polygons.resize( _polygons.size() + 1 );
3252         _polygons[ iPolygon ]._polyLinks.reserve( 20 );
3253         _polygons[ iPolygon ]._links.reserve( 20 );
3254       }
3255       _Face& polygon = _polygons[ iPolygon ];
3256
3257       _OrientedLink* curLink = 0;
3258       _Node*         curNode;
3259       if (( !hasEdgeIntersections ) ||
3260           ( nbFreeLinks < 4 && nbVertexNodes == 0 ))
3261       {
3262         // get a remaining link to start from
3263         for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
3264           if (( curLink = freeLinks[ iL ] ))
3265             freeLinks[ iL ] = 0;
3266         polygon._links.push_back( *curLink );
3267         --nbFreeLinks;
3268         do
3269         {
3270           // find all links connected to curLink
3271           curNode = curLink->FirstNode();
3272           curLink = 0;
3273           for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
3274             if ( freeLinks[ iL ] && freeLinks[ iL ]->LastNode() == curNode )
3275             {
3276               curLink = freeLinks[ iL ];
3277               freeLinks[ iL ] = 0;
3278               --nbFreeLinks;
3279               polygon._links.push_back( *curLink );
3280             }
3281         } while ( curLink );
3282       }
3283       else // there are intersections with EDGEs
3284       {
3285         // get a remaining link to start from, one lying on minimal nb of FACEs
3286         {
3287           typedef pair< TGeomID, int > TFaceOfLink;
3288           TFaceOfLink faceOfLink( -1, -1 );
3289           TFaceOfLink facesOfLink[3] = { faceOfLink, faceOfLink, faceOfLink };
3290           for ( size_t iL = 0; iL < freeLinks.size(); ++iL )
3291             if ( freeLinks[ iL ] )
3292             {
3293               faces = freeLinks[ iL ]->GetNotUsedFace( usedFaceIDs );
3294               if ( faces.size() == 1 )
3295               {
3296                 faceOfLink = TFaceOfLink( faces[0], iL );
3297                 if ( !freeLinks[ iL ]->HasEdgeNodes() )
3298                   break;
3299                 facesOfLink[0] = faceOfLink;
3300               }
3301               else if ( facesOfLink[0].first < 0 )
3302               {
3303                 faceOfLink = TFaceOfLink(( faces.empty() ? -1 : faces[0]), iL );
3304                 facesOfLink[ 1 + faces.empty() ] = faceOfLink;
3305               }
3306             }
3307           for ( int i = 0; faceOfLink.first < 0 && i < 3; ++i )
3308             faceOfLink = facesOfLink[i];
3309
3310           if ( faceOfLink.first < 0 ) // all faces used
3311           {
3312             for ( size_t iL = 0; iL < freeLinks.size() && faceOfLink.first < 1; ++iL )
3313               if (( curLink = freeLinks[ iL ]))
3314               {
3315                 faceOfLink.first = 
3316                   curLink->FirstNode()->IsLinked( curLink->LastNode()->_intPoint );
3317                 faceOfLink.second = iL;
3318               }
3319             usedFaceIDs.clear();
3320           }
3321           curFace = faceOfLink.first;
3322           curLink = freeLinks[ faceOfLink.second ];
3323           freeLinks[ faceOfLink.second ] = 0;
3324         }
3325         usedFaceIDs.insert( curFace );
3326         polygon._links.push_back( *curLink );
3327         --nbFreeLinks;
3328
3329         // find all links lying on a curFace
3330         do
3331         {
3332           // go forward from curLink
3333           curNode = curLink->LastNode();
3334           curLink = 0;
3335           for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
3336             if ( freeLinks[ iL ] &&
3337                  freeLinks[ iL ]->FirstNode() == curNode &&
3338                  freeLinks[ iL ]->LastNode()->IsOnFace( curFace ))
3339             {
3340               curLink = freeLinks[ iL ];
3341               freeLinks[ iL ] = 0;
3342               polygon._links.push_back( *curLink );
3343               --nbFreeLinks;
3344             }
3345         } while ( curLink );
3346
3347         std::reverse( polygon._links.begin(), polygon._links.end() );
3348
3349         curLink = & polygon._links.back();
3350         do
3351         {
3352           // go backward from curLink
3353           curNode = curLink->FirstNode();
3354           curLink = 0;
3355           for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
3356             if ( freeLinks[ iL ] &&
3357                  freeLinks[ iL ]->LastNode() == curNode &&
3358                  freeLinks[ iL ]->FirstNode()->IsOnFace( curFace ))
3359             {
3360               curLink = freeLinks[ iL ];
3361               freeLinks[ iL ] = 0;
3362               polygon._links.push_back( *curLink );
3363               --nbFreeLinks;
3364             }
3365         } while ( curLink );
3366
3367         curNode = polygon._links.back().FirstNode();
3368
3369         if ( polygon._links[0].LastNode() != curNode )
3370         {
3371           if ( nbVertexNodes > 0 )
3372           {
3373             // add links with _vIntNodes if not already used
3374             chainNodes.clear();
3375             for ( size_t iN = 0; iN < _vIntNodes.size(); ++iN )
3376               if ( !_vIntNodes[ iN ]->IsUsedInFace() &&
3377                    _vIntNodes[ iN ]->IsOnFace( curFace ))
3378               {
3379                 _vIntNodes[ iN ]->_usedInFace = &polygon;
3380                 chainNodes.push_back( _vIntNodes[ iN ] );
3381               }
3382             if ( chainNodes.size() > 1 &&
3383                  curFace != _grid->PseudoIntExtFaceID() ) /////// TODO
3384             {
3385               sortVertexNodes( chainNodes, curNode, curFace );
3386             }
3387             for ( size_t i = 0; i < chainNodes.size(); ++i )
3388             {
3389               polygon.AddPolyLink( chainNodes[ i ], curNode );
3390               curNode = chainNodes[ i ];
3391               freeLinks.push_back( &polygon._links.back() );
3392               ++nbFreeLinks;
3393             }
3394             nbVertexNodes -= chainNodes.size();
3395           }
3396           // if ( polygon._links.size() > 1 )
3397           {
3398             polygon.AddPolyLink( polygon._links[0].LastNode(), curNode );
3399             freeLinks.push_back( &polygon._links.back() );
3400             ++nbFreeLinks;
3401           }
3402         }
3403       } // if there are intersections with EDGEs
3404
3405       if ( polygon._links.size() < 2 ||
3406            polygon._links[0].LastNode() != polygon._links.back().FirstNode() )
3407       {
3408         _polygons.clear();
3409         break; // closed polygon not found -> invalid polyhedron
3410       }
3411
3412       if ( polygon._links.size() == 2 )
3413       {
3414         if ( freeLinks.back() == &polygon._links.back() )
3415         {
3416           freeLinks.pop_back();
3417           --nbFreeLinks;
3418         }
3419         if ( polygon._links.front().NbFaces() > 0 )
3420           polygon._links.back().AddFace( polygon._links.front()._link->_faces[0] );
3421         if ( polygon._links.back().NbFaces() > 0 )
3422           polygon._links.front().AddFace( polygon._links.back()._link->_faces[0] );
3423
3424         if ( iPolygon == _polygons.size()-1 )
3425           _polygons.pop_back();
3426       }
3427       else // polygon._links.size() >= 2
3428       {
3429         // add polygon to its links
3430         for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
3431         {
3432           polygon._links[ iL ].AddFace( &polygon );
3433           polygon._links[ iL ].Reverse();
3434         }
3435         if ( /*hasEdgeIntersections &&*/ iPolygon == _polygons.size() - 1 )
3436         {
3437           // check that a polygon does not lie on a hexa side
3438           coplanarPolyg = 0;
3439           for ( size_t iL = 0; iL < polygon._links.size() && !coplanarPolyg; ++iL )
3440           {
3441             if ( polygon._links[ iL ].NbFaces() < 2 )
3442               continue; // it's a just added free link
3443             // look for a polygon made on a hexa side and sharing
3444             // two or more haxa links
3445             size_t iL2;
3446             coplanarPolyg = polygon._links[ iL ]._link->_faces[0];
3447             for ( iL2 = iL + 1; iL2 < polygon._links.size(); ++iL2 )
3448               if ( polygon._links[ iL2 ]._link->_faces[0] == coplanarPolyg &&
3449                    !coplanarPolyg->IsPolyLink( polygon._links[ iL  ]) &&
3450                    !coplanarPolyg->IsPolyLink( polygon._links[ iL2 ]) &&
3451                    coplanarPolyg < & _polygons[ nbQuadPolygons ])
3452                 break;
3453             if ( iL2 == polygon._links.size() )
3454               coplanarPolyg = 0;
3455           }
3456           if ( coplanarPolyg ) // coplanar polygon found
3457           {
3458             freeLinks.resize( freeLinks.size() - polygon._polyLinks.size() );
3459             nbFreeLinks -= polygon._polyLinks.size();
3460
3461             // an E_IntersectPoint used to mark nodes of coplanarPolyg
3462             // as lying on curFace while they are not at intersection with geometry
3463             ipTmp._faceIDs.resize(1);
3464             ipTmp._faceIDs[0] = curFace;
3465
3466             // fill freeLinks with links not shared by coplanarPolyg and polygon
3467             for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
3468               if ( polygon._links[ iL ]._link->_faces[1] &&
3469                    polygon._links[ iL ]._link->_faces[0] != coplanarPolyg )
3470               {
3471                 _Face* p = polygon._links[ iL ]._link->_faces[0];
3472                 for ( size_t iL2 = 0; iL2 < p->_links.size(); ++iL2 )
3473                   if ( p->_links[ iL2 ]._link == polygon._links[ iL ]._link )
3474                   {
3475                     freeLinks.push_back( & p->_links[ iL2 ] );
3476                     ++nbFreeLinks;
3477                     freeLinks.back()->RemoveFace( &polygon );
3478                     break;
3479                   }
3480               }
3481             for ( size_t iL = 0; iL < coplanarPolyg->_links.size(); ++iL )
3482               if ( coplanarPolyg->_links[ iL ]._link->_faces[1] &&
3483                    coplanarPolyg->_links[ iL ]._link->_faces[1] != &polygon )
3484               {
3485                 _Face* p = coplanarPolyg->_links[ iL ]._link->_faces[0];
3486                 if ( p == coplanarPolyg )
3487                   p = coplanarPolyg->_links[ iL ]._link->_faces[1];
3488                 for ( size_t iL2 = 0; iL2 < p->_links.size(); ++iL2 )
3489                   if ( p->_links[ iL2 ]._link == coplanarPolyg->_links[ iL ]._link )
3490                   {
3491                     // set links of coplanarPolyg in place of used freeLinks
3492                     // to re-create coplanarPolyg next
3493                     size_t iL3 = 0;
3494                     for ( ; iL3 < freeLinks.size() && freeLinks[ iL3 ]; ++iL3 );
3495                     if ( iL3 < freeLinks.size() )
3496                       freeLinks[ iL3 ] = ( & p->_links[ iL2 ] );
3497                     else
3498                       freeLinks.push_back( & p->_links[ iL2 ] );
3499                     ++nbFreeLinks;
3500                     freeLinks[ iL3 ]->RemoveFace( coplanarPolyg );
3501                     //  mark nodes of coplanarPolyg as lying on curFace
3502                     for ( int iN = 0; iN < 2; ++iN )
3503                     {
3504                       _Node* n = freeLinks[ iL3 ]->_link->_nodes[ iN ];
3505                       bool added = false;
3506                       if ( n->_intPoint ) added = n->_intPoint->Add( ipTmp._faceIDs );
3507                       else                        n->_intPoint = &ipTmp;
3508                       if ( added )
3509                         tmpAddedFace[ ipTmp._faceIDs[0] ].push_back( n->_intPoint );
3510                     }
3511                     break;
3512                   }
3513               }
3514             // set coplanarPolyg to be re-created next
3515             for ( size_t iP = 0; iP < _polygons.size(); ++iP )
3516               if ( coplanarPolyg == & _polygons[ iP ] )
3517               {
3518                 iPolygon = iP;
3519                 _polygons[ iPolygon ]._links.clear();
3520                 _polygons[ iPolygon ]._polyLinks.clear();
3521                 break;
3522               }
3523             _polygons.pop_back();
3524             usedFaceIDs.erase( curFace );
3525             continue;
3526           } // if ( coplanarPolyg )
3527         } // if ( hasEdgeIntersections ) - search for coplanarPolyg
3528
3529         iPolygon = _polygons.size();
3530
3531       } // end of case ( polygon._links.size() > 2 )
3532     } // while ( nbFreeLinks > 0 )
3533
3534     for ( auto & face_ip : tmpAddedFace )
3535     {
3536       curFace = face_ip.first;
3537       for ( const B_IntersectPoint* ip : face_ip.second )
3538       {
3539         auto it = std::find( ip->_faceIDs.begin(), ip->_faceIDs.end(), curFace );
3540         if ( it != ip->_faceIDs.end() )
3541           ip->_faceIDs.erase( it );
3542       }
3543     }
3544
3545     if ( _polygons.size() < 3 )
3546       return false;
3547
3548     // check volume size
3549     double volSize = 0;
3550     _hasTooSmall = ! checkPolyhedronSize( intFlag & IS_CUT_BY_INTERNAL_FACE, volSize );
3551
3552     for ( size_t i = 0; i < 8; ++i )
3553       if ( _hexNodes[ i ]._intPoint == &ipTmp )
3554         _hexNodes[ i ]._intPoint = 0;
3555
3556     if ( _hasTooSmall )
3557       return false; // too small volume
3558
3559
3560     // Try to find out names of no-name polygons (issue # 19887)
3561     if ( _grid->IsToRemoveExcessEntities() && _polygons.back()._name == SMESH_Block::ID_NONE )
3562     {
3563       gp_XYZ uvwCenter =
3564         0.5 * ( _grid->_coords[0][_i] + _grid->_coords[0][_i+1] ) * _grid->_axes[0] +
3565         0.5 * ( _grid->_coords[1][_j] + _grid->_coords[1][_j+1] ) * _grid->_axes[1] +
3566         0.5 * ( _grid->_coords[2][_k] + _grid->_coords[2][_k+1] ) * _grid->_axes[2];
3567       for ( size_t i = _polygons.size() - 1; _polygons[i]._name == SMESH_Block::ID_NONE; --i )
3568       {
3569         _Face& face = _polygons[ i ];
3570         Bnd_Box bb;
3571         gp_Pnt uvw;
3572         for ( size_t iL = 0; iL < face._links.size(); ++iL )
3573         {
3574           _Node* n = face._links[ iL ].FirstNode();
3575           gp_XYZ p = SMESH_NodeXYZ( n->Node() );
3576           _grid->ComputeUVW( p, uvw.ChangeCoord().ChangeData() );
3577           bb.Add( uvw );
3578         }
3579         gp_Pnt pMin = bb.CornerMin();
3580         if ( bb.IsXThin( _grid->_tol ))
3581           face._name = pMin.X() < uvwCenter.X() ? SMESH_Block::ID_F0yz : SMESH_Block::ID_F1yz;
3582         else if ( bb.IsYThin( _grid->_tol ))
3583           face._name = pMin.Y() < uvwCenter.Y() ? SMESH_Block::ID_Fx0z : SMESH_Block::ID_Fx1z;
3584         else if ( bb.IsZThin( _grid->_tol ))
3585           face._name = pMin.Z() < uvwCenter.Z() ? SMESH_Block::ID_Fxy0 : SMESH_Block::ID_Fxy1;
3586       }
3587     }
3588
3589     _volumeDefs._nodes.clear();
3590     _volumeDefs._quantities.clear();
3591     _volumeDefs._names.clear();
3592     // create a classic cell if possible
3593
3594     int nbPolygons = 0;
3595     for ( size_t iF = 0; iF < _polygons.size(); ++iF )
3596       nbPolygons += (_polygons[ iF ]._links.size() > 2 );
3597
3598     //const int nbNodes = _nbCornerNodes + nbIntersections;
3599     int nbNodes = 0;
3600     for ( size_t i = 0; i < 8; ++i )
3601       nbNodes += _hexNodes[ i ].IsUsedInFace();
3602     for ( size_t i = 0; i < _intNodes.size(); ++i )
3603       nbNodes += _intNodes[ i ].IsUsedInFace();
3604
3605     bool isClassicElem = false;
3606     if (      nbNodes == 8 && nbPolygons == 6 ) isClassicElem = addHexa();
3607     else if ( nbNodes == 4 && nbPolygons == 4 ) isClassicElem = addTetra();
3608     else if ( nbNodes == 6 && nbPolygons == 5 ) isClassicElem = addPenta();
3609     else if ( nbNodes == 5 && nbPolygons == 5 ) isClassicElem = addPyra ();
3610     if ( !isClassicElem )
3611     {
3612       for ( size_t iF = 0; iF < _polygons.size(); ++iF )
3613       {
3614         const size_t nbLinks = _polygons[ iF ]._links.size();
3615         if ( nbLinks < 3 ) continue;
3616         _volumeDefs._quantities.push_back( nbLinks );
3617         _volumeDefs._names.push_back( _polygons[ iF ]._name );
3618         for ( size_t iL = 0; iL < nbLinks; ++iL )
3619           _volumeDefs._nodes.push_back( _polygons[ iF ]._links[ iL ].FirstNode() );
3620       }
3621     }
3622     _volumeDefs._solidID = solid->ID();
3623     _volumeDefs._size    = volSize;
3624
3625     return !_volumeDefs._nodes.empty();
3626   }
3627
3628   template<typename Type>
3629   void computeHexa(Type& hex)
3630   {
3631     if ( hex ) 
3632       hex->computeElements();
3633   }
3634
3635   // Implement parallel computation of Hexa with c++ thread implementation
3636   template<typename Iterator, class Function>
3637   void parallel_for(const Iterator& first, const Iterator& last, Function&& f, const int nthreads = 1)
3638   {
3639       const unsigned int group = ((last-first))/std::abs(nthreads);
3640
3641       std::vector<std::thread> threads;
3642       threads.reserve(nthreads);
3643       Iterator it = first;
3644       for (; it < last-group; it += group) {
3645           // to create a thread 
3646           // Pass iterators by value and the function by reference!
3647           auto lambda = [=,&f](){ std::for_each(it, std::min(it+group, last), f);};
3648
3649           // stack the threads 
3650           threads.push_back( std::thread( lambda ) );
3651       }
3652
3653       std::for_each(it, last, f); // last steps while we wait for other threads
3654       std::for_each(threads.begin(), threads.end(), [](std::thread& x){x.join();});
3655   }
3656   //================================================================================
3657   /*!
3658    * \brief Create elements in the mesh
3659    */
3660   int Hexahedron::MakeElements(SMESH_MesherHelper&                      helper,
3661                                const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap)
3662   {
3663     SMESHDS_Mesh* mesh = helper.GetMeshDS();
3664
3665     CellsAroundLink c( _grid, 0 );
3666     const size_t nbGridCells = c._nbCells[0] * c._nbCells[1] * c._nbCells[2];
3667     vector< Hexahedron* > allHexa( nbGridCells, 0 );
3668     int nbIntHex = 0;
3669
3670     // set intersection nodes from GridLine's to links of allHexa
3671     int i,j,k, cellIndex, iLink;
3672     for ( int iDir = 0; iDir < 3; ++iDir )
3673     {
3674       // loop on GridLine's parallel to iDir
3675       LineIndexer lineInd = _grid->GetLineIndexer( iDir );
3676       CellsAroundLink fourCells( _grid, iDir );
3677       for ( ; lineInd.More(); ++lineInd )
3678       {
3679         GridLine& line = _grid->_lines[ iDir ][ lineInd.LineIndex() ];
3680         multiset< F_IntersectPoint >::const_iterator ip = line._intPoints.begin();
3681         for ( ; ip != line._intPoints.end(); ++ip )
3682         {
3683           // if ( !ip->_node ) continue; // intersection at a grid node
3684           lineInd.SetIndexOnLine( ip->_indexOnLine );
3685           fourCells.Init( lineInd.I(), lineInd.J(), lineInd.K() );
3686           for ( int iL = 0; iL < 4; ++iL ) // loop on 4 cells sharing a link
3687           {
3688             if ( !fourCells.GetCell( iL, i,j,k, cellIndex, iLink ))
3689               continue;
3690             Hexahedron *& hex = allHexa[ cellIndex ];
3691             if ( !hex)
3692             {
3693               hex = new Hexahedron( *this, i, j, k, cellIndex );
3694               ++nbIntHex;
3695             }
3696             hex->_hexLinks[iLink]._fIntPoints.push_back( &(*ip) );
3697             hex->_nbFaceIntNodes += bool( ip->_node );
3698           }
3699         }
3700       }
3701     }
3702
3703     // implement geom edges into the mesh
3704     addEdges( helper, allHexa, edge2faceIDsMap );
3705
3706     // add not split hexahedra to the mesh
3707     int nbAdded = 0;
3708     TGeomID solidIDs[20];
3709     vector< Hexahedron* > intHexa; intHexa.reserve( nbIntHex );
3710     vector< const SMDS_MeshElement* > boundaryVolumes; boundaryVolumes.reserve( nbIntHex * 1.1 );
3711     for ( size_t i = 0; i < allHexa.size(); ++i )
3712     {
3713       // initialize this by not cut allHexa[ i ]
3714       Hexahedron * & hex = allHexa[ i ];
3715       if ( hex ) // split hexahedron
3716       {
3717         intHexa.push_back( hex );
3718         if ( hex->_nbFaceIntNodes > 0 ||
3719              hex->_eIntPoints.size() > 0 ||
3720              hex->getSolids( solidIDs ) > 1 )
3721           continue; // treat intersected hex later in parallel
3722         this->init( hex->_i, hex->_j, hex->_k );
3723       }
3724       else
3725       {
3726         this->init( i ); // == init(i,j,k)
3727       }
3728       if (( _nbCornerNodes == 8 ) &&
3729           ( _nbBndNodes < _nbCornerNodes || !isInHole() ))
3730       {
3731         // order of _hexNodes is defined by enum SMESH_Block::TShapeID
3732         SMDS_MeshElement* el =
3733           mesh->AddVolume( _hexNodes[0].Node(), _hexNodes[2].Node(),
3734                            _hexNodes[3].Node(), _hexNodes[1].Node(),
3735                            _hexNodes[4].Node(), _hexNodes[6].Node(),
3736                            _hexNodes[7].Node(), _hexNodes[5].Node() );
3737         TGeomID solidID = 0;
3738         if ( _nbBndNodes < _nbCornerNodes )
3739         {
3740           for ( int iN = 0; iN < 8 &&  !solidID; ++iN )
3741             if ( !_hexNodes[iN]._intPoint ) // no intersection
3742               solidID = _hexNodes[iN].Node()->GetShapeID();
3743         }
3744         else
3745         {
3746           getSolids( solidIDs );
3747           solidID = solidIDs[0];
3748         }
3749         mesh->SetMeshElementOnShape( el, solidID );
3750         ++nbAdded;
3751         if ( hex )
3752           intHexa.pop_back();
3753         if ( _grid->_toCreateFaces && _nbBndNodes >= 3 )
3754         {
3755           boundaryVolumes.push_back( el );
3756           el->setIsMarked( true );
3757         }
3758       }
3759       else if ( _nbCornerNodes > 3 && !hex )
3760       {
3761         // all intersections of hex with geometry are at grid nodes
3762         hex = new Hexahedron( *this, _i, _j, _k, i );
3763         intHexa.push_back( hex );
3764       }
3765     }
3766
3767     // compute definitions of volumes resulted from hexadron intersection
3768 #ifdef WITH_TBB
3769     auto numOfThreads = std::thread::hardware_concurrency();
3770     numOfThreads = (numOfThreads != 0) ? numOfThreads : 1;
3771     parallel_for(intHexa.begin(), intHexa.end(), computeHexa<Hexahedron*>, numOfThreads ); 
3772 #else
3773     for ( size_t i = 0; i < intHexa.size(); ++i )
3774       if ( Hexahedron * hex = intHexa[ i ] )
3775         hex->computeElements();
3776 #endif
3777
3778     // simplify polyhedrons
3779     if ( _grid->IsToRemoveExcessEntities() )
3780     {
3781       for ( size_t i = 0; i < intHexa.size(); ++i )
3782         if ( Hexahedron * hex = intHexa[ i ] )
3783           hex->removeExcessSideDivision( allHexa );
3784
3785       for ( size_t i = 0; i < intHexa.size(); ++i )
3786         if ( Hexahedron * hex = intHexa[ i ] )
3787           hex->removeExcessNodes( allHexa );
3788     }
3789
3790     // add volumes
3791     for ( size_t i = 0; i < intHexa.size(); ++i )
3792       if ( Hexahedron * hex = intHexa[ i ] )
3793         nbAdded += hex->addVolumes( helper );
3794
3795     // fill boundaryVolumes with volumes neighboring too small skipped volumes
3796     if ( _grid->_toCreateFaces )
3797     {
3798       for ( size_t i = 0; i < intHexa.size(); ++i )
3799         if ( Hexahedron * hex = intHexa[ i ] )
3800           hex->getBoundaryElems( boundaryVolumes );
3801     }
3802
3803     // merge nodes on outer sub-shapes with pre-existing ones
3804     TopTools_DataMapIteratorOfDataMapOfShapeInteger s2nIt( _grid->_geometry._shape2NbNodes );
3805     for ( ; s2nIt.More(); s2nIt.Next() )
3806       if ( s2nIt.Value() > 0 )
3807         if ( SMESHDS_SubMesh* sm = mesh->MeshElements( s2nIt.Key() ))
3808         {
3809           TIDSortedNodeSet smNodes( SMDS_MeshElement::iterator( sm->GetNodes() ),
3810                                     SMDS_MeshElement::iterator() );
3811           SMESH_MeshEditor::TListOfListOfNodes equalNodes;
3812           SMESH_MeshEditor editor( helper.GetMesh() );
3813           editor.FindCoincidentNodes( smNodes, 10 * _grid->_tol, equalNodes,
3814                                       /*SeparateCornersAndMedium =*/ false);
3815           if ((int) equalNodes.size() <= s2nIt.Value() )
3816             editor.MergeNodes( equalNodes );
3817         }
3818
3819     // create boundary mesh faces
3820     addFaces( helper, boundaryVolumes );
3821
3822     // create mesh edges
3823     addSegments( helper, edge2faceIDsMap );
3824
3825     for ( size_t i = 0; i < allHexa.size(); ++i )
3826       if ( allHexa[ i ] )
3827         delete allHexa[ i ];
3828
3829     return nbAdded;
3830   }
3831
3832   //================================================================================
3833   /*!
3834    * \brief Implements geom edges into the mesh
3835    */
3836   void Hexahedron::addEdges(SMESH_MesherHelper&                      helper,
3837                             vector< Hexahedron* >&                   hexes,
3838                             const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap)
3839   {
3840     if ( edge2faceIDsMap.empty() ) return;
3841
3842     // Prepare planes for intersecting with EDGEs
3843     GridPlanes pln[3];
3844     {
3845       for ( int iDirZ = 0; iDirZ < 3; ++iDirZ ) // iDirZ gives normal direction to planes
3846       {
3847         GridPlanes& planes = pln[ iDirZ ];
3848         int iDirX = ( iDirZ + 1 ) % 3;
3849         int iDirY = ( iDirZ + 2 ) % 3;
3850         planes._zNorm  = ( _grid->_axes[ iDirX ] ^ _grid->_axes[ iDirY ] ).Normalized();
3851         planes._zProjs.resize ( _grid->_coords[ iDirZ ].size() );
3852         planes._zProjs [0] = 0;
3853         const double       zFactor = _grid->_axes[ iDirZ ] * planes._zNorm;
3854         const vector< double > & u = _grid->_coords[ iDirZ ];
3855         for ( size_t i = 1; i < planes._zProjs.size(); ++i )
3856         {
3857           planes._zProjs [i] = zFactor * ( u[i] - u[0] );
3858         }
3859       }
3860     }
3861     const double deflection = _grid->_minCellSize / 20.;
3862     const double tol        = _grid->_tol;
3863     E_IntersectPoint ip;
3864
3865     TColStd_MapOfInteger intEdgeIDs; // IDs of not shared INTERNAL EDGES
3866
3867     // Intersect EDGEs with the planes
3868     map< TGeomID, vector< TGeomID > >::const_iterator e2fIt = edge2faceIDsMap.begin();
3869     for ( ; e2fIt != edge2faceIDsMap.end(); ++e2fIt )
3870     {
3871       const TGeomID  edgeID = e2fIt->first;
3872       const TopoDS_Edge & E = TopoDS::Edge( _grid->Shape( edgeID ));
3873       BRepAdaptor_Curve curve( E );
3874       TopoDS_Vertex v1 = helper.IthVertex( 0, E, false );
3875       TopoDS_Vertex v2 = helper.IthVertex( 1, E, false );
3876
3877       ip._faceIDs = e2fIt->second;
3878       ip._shapeID = edgeID;
3879
3880       bool isInternal = ( ip._faceIDs.size() == 1 && _grid->IsInternal( edgeID ));
3881       if ( isInternal )
3882       {
3883         intEdgeIDs.Add( edgeID );
3884         intEdgeIDs.Add( _grid->ShapeID( v1 ));
3885         intEdgeIDs.Add( _grid->ShapeID( v2 ));
3886       }
3887
3888       // discretize the EDGE
3889       GCPnts_UniformDeflection discret( curve, deflection, true );
3890       if ( !discret.IsDone() || discret.NbPoints() < 2 )
3891         continue;
3892
3893       // perform intersection
3894       E_IntersectPoint* eip, *vip = 0;
3895       for ( int iDirZ = 0; iDirZ < 3; ++iDirZ )
3896       {
3897         GridPlanes& planes = pln[ iDirZ ];
3898         int      iDirX = ( iDirZ + 1 ) % 3;
3899         int      iDirY = ( iDirZ + 2 ) % 3;
3900         double    xLen = _grid->_coords[ iDirX ].back() - _grid->_coords[ iDirX ][0];
3901         double    yLen = _grid->_coords[ iDirY ].back() - _grid->_coords[ iDirY ][0];
3902         double    zLen = _grid->_coords[ iDirZ ].back() - _grid->_coords[ iDirZ ][0];
3903         int dIJK[3], d000[3] = { 0,0,0 };
3904         double o[3] = { _grid->_coords[0][0],
3905                         _grid->_coords[1][0],
3906                         _grid->_coords[2][0] };
3907
3908         // locate the 1st point of a segment within the grid
3909         gp_XYZ p1     = discret.Value( 1 ).XYZ();
3910         double u1     = discret.Parameter( 1 );
3911         double zProj1 = planes._zNorm * ( p1 - _grid->_origin );
3912
3913         _grid->ComputeUVW( p1, ip._uvw );
3914         int iX1 = int(( ip._uvw[iDirX] - o[iDirX]) / xLen * (_grid->_coords[ iDirX ].size() - 1));
3915         int iY1 = int(( ip._uvw[iDirY] - o[iDirY]) / yLen * (_grid->_coords[ iDirY ].size() - 1));
3916         int iZ1 = int(( ip._uvw[iDirZ] - o[iDirZ]) / zLen * (_grid->_coords[ iDirZ ].size() - 1));
3917         locateValue( iX1, ip._uvw[iDirX], _grid->_coords[ iDirX ], dIJK[ iDirX ], tol );
3918         locateValue( iY1, ip._uvw[iDirY], _grid->_coords[ iDirY ], dIJK[ iDirY ], tol );
3919         locateValue( iZ1, ip._uvw[iDirZ], _grid->_coords[ iDirZ ], dIJK[ iDirZ ], tol );
3920
3921         int ijk[3]; // grid index where a segment intersects a plane
3922         ijk[ iDirX ] = iX1;
3923         ijk[ iDirY ] = iY1;
3924         ijk[ iDirZ ] = iZ1;
3925
3926         // add the 1st vertex point to a hexahedron
3927         if ( iDirZ == 0 )
3928         {
3929           ip._point   = p1;
3930           ip._shapeID = _grid->ShapeID( v1 );
3931           vip = _grid->Add( ip );
3932           _grid->UpdateFacesOfVertex( *vip, v1 );
3933           if ( isInternal )
3934             vip->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
3935           if ( !addIntersection( vip, hexes, ijk, d000 ))
3936             _grid->Remove( vip );
3937           ip._shapeID = edgeID;
3938         }
3939         for ( int iP = 2; iP <= discret.NbPoints(); ++iP )
3940         {
3941           // locate the 2nd point of a segment within the grid
3942           gp_XYZ p2     = discret.Value( iP ).XYZ();
3943           double u2     = discret.Parameter( iP );
3944           double zProj2 = planes._zNorm * ( p2 - _grid->_origin );
3945           int    iZ2    = iZ1;
3946           if ( Abs( zProj2 - zProj1 ) > std::numeric_limits<double>::min() )
3947           {
3948             locateValue( iZ2, zProj2, planes._zProjs, dIJK[ iDirZ ], tol );
3949
3950             // treat intersections with planes between 2 end points of a segment
3951             int dZ = ( iZ1 <= iZ2 ) ? +1 : -1;
3952             int iZ = iZ1 + ( iZ1 < iZ2 );
3953             for ( int i = 0, nb = Abs( iZ1 - iZ2 ); i < nb; ++i, iZ += dZ )
3954             {
3955               ip._point = findIntPoint( u1, zProj1, u2, zProj2,
3956                                         planes._zProjs[ iZ ],
3957                                         curve, planes._zNorm, _grid->_origin );
3958               _grid->ComputeUVW( ip._point.XYZ(), ip._uvw );
3959               locateValue( ijk[iDirX], ip._uvw[iDirX], _grid->_coords[iDirX], dIJK[iDirX], tol );
3960               locateValue( ijk[iDirY], ip._uvw[iDirY], _grid->_coords[iDirY], dIJK[iDirY], tol );
3961               ijk[ iDirZ ] = iZ;
3962
3963               // add ip to hex "above" the plane
3964               eip = _grid->Add( ip );
3965               if ( isInternal )
3966                 eip->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
3967               dIJK[ iDirZ ] = 0;
3968               bool added = addIntersection( eip, hexes, ijk, dIJK);
3969
3970               // add ip to hex "below" the plane
3971               ijk[ iDirZ ] = iZ-1;
3972               if ( !addIntersection( eip, hexes, ijk, dIJK ) &&
3973                    !added )
3974                 _grid->Remove( eip );
3975             }
3976           }
3977           iZ1    = iZ2;
3978           p1     = p2;
3979           u1     = u2;
3980           zProj1 = zProj2;
3981         }
3982         // add the 2nd vertex point to a hexahedron
3983         if ( iDirZ == 0 )
3984         {
3985           ip._point   = p1;
3986           ip._shapeID = _grid->ShapeID( v2 );
3987           _grid->ComputeUVW( p1, ip._uvw );
3988           locateValue( ijk[iDirX], ip._uvw[iDirX], _grid->_coords[iDirX], dIJK[iDirX], tol );
3989           locateValue( ijk[iDirY], ip._uvw[iDirY], _grid->_coords[iDirY], dIJK[iDirY], tol );
3990           ijk[ iDirZ ] = iZ1;
3991           bool sameV = ( v1.IsSame( v2 ));
3992           if ( !sameV )
3993           {
3994             vip = _grid->Add( ip );
3995             _grid->UpdateFacesOfVertex( *vip, v2 );
3996             if ( isInternal )
3997               vip->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
3998           }
3999           if ( !addIntersection( vip, hexes, ijk, d000 ) && !sameV )
4000             _grid->Remove( vip );
4001           ip._shapeID = edgeID;
4002         }
4003       } // loop on 3 grid directions
4004     } // loop on EDGEs
4005
4006
4007     if ( intEdgeIDs.Size() > 0 )
4008       cutByExtendedInternal( hexes, intEdgeIDs );
4009
4010     return;
4011   }
4012
4013   //================================================================================
4014   /*!
4015    * \brief Fully cut hexes that are partially cut by INTERNAL FACE.
4016    *        Cut them by extended INTERNAL FACE.
4017    */
4018   void Hexahedron::cutByExtendedInternal( std::vector< Hexahedron* >& hexes,
4019                                           const TColStd_MapOfInteger& intEdgeIDs )
4020   {
4021     IntAna_IntConicQuad intersection;
4022     SMESHDS_Mesh* meshDS = _grid->_helper->GetMeshDS();
4023     const double tol2 = _grid->_tol * _grid->_tol;
4024
4025     for ( size_t iH = 0; iH < hexes.size(); ++iH )
4026     {
4027       Hexahedron* hex = hexes[ iH ];
4028       if ( !hex || hex->_eIntPoints.size() < 2 )
4029         continue;
4030       if ( !intEdgeIDs.Contains( hex->_eIntPoints.back()->_shapeID ))
4031         continue;
4032
4033       // get 3 points on INTERNAL FACE to construct a cutting plane
4034       gp_Pnt p1 = hex->_eIntPoints[0]->_point;
4035       gp_Pnt p2 = hex->_eIntPoints[1]->_point;
4036       gp_Pnt p3 = hex->mostDistantInternalPnt( iH, p1, p2 );
4037
4038       gp_Vec norm = gp_Vec( p1, p2 ) ^ gp_Vec( p1, p3 );
4039       gp_Pln pln;
4040       try {
4041         pln = gp_Pln( p1, norm );
4042       }
4043       catch(...)
4044       {
4045         continue;
4046       }
4047
4048       TGeomID intFaceID = hex->_eIntPoints.back()->_faceIDs.front(); // FACE being "extended"
4049       TGeomID   solidID = _grid->GetSolid( intFaceID )->ID();
4050
4051       // cut links by the plane
4052       //bool isCut = false;
4053       for ( int iLink = 0; iLink < 12; ++iLink )
4054       {
4055         _Link& link = hex->_hexLinks[ iLink ];
4056         if ( !link._fIntPoints.empty() )
4057         {
4058           // if ( link._fIntPoints[0]->_faceIDs.back() == _grid->PseudoIntExtFaceID() )
4059           //   isCut = true;
4060           continue; // already cut link
4061         }
4062         if ( !link._nodes[0]->Node() ||
4063              !link._nodes[1]->Node() )
4064           continue; // outside link
4065
4066         if ( link._nodes[0]->IsOnFace( intFaceID ))
4067         {
4068           if ( link._nodes[0]->_intPoint->_faceIDs.back() != _grid->PseudoIntExtFaceID() )
4069             if ( p1.SquareDistance( link._nodes[0]->Point() ) < tol2  ||
4070                  p2.SquareDistance( link._nodes[0]->Point() ) < tol2 )
4071               link._nodes[0]->_intPoint->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
4072           continue; // link is cut by FACE being "extended"
4073         }
4074         if ( link._nodes[1]->IsOnFace( intFaceID ))
4075         {
4076           if ( link._nodes[1]->_intPoint->_faceIDs.back() != _grid->PseudoIntExtFaceID() )
4077             if ( p1.SquareDistance( link._nodes[1]->Point() ) < tol2  ||
4078                  p2.SquareDistance( link._nodes[1]->Point() ) < tol2 )
4079               link._nodes[1]->_intPoint->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
4080           continue; // link is cut by FACE being "extended"
4081         }
4082         gp_Pnt p4 = link._nodes[0]->Point();
4083         gp_Pnt p5 = link._nodes[1]->Point();
4084         gp_Lin line( p4, gp_Vec( p4, p5 ));
4085
4086         intersection.Perform( line, pln );
4087         if ( !intersection.IsDone() ||
4088              intersection.IsInQuadric() ||
4089              intersection.IsParallel() ||
4090              intersection.NbPoints() < 1 )
4091           continue;
4092
4093         double u = intersection.ParamOnConic(1);
4094         if ( u + _grid->_tol < 0 )
4095           continue;
4096         int       iDir = iLink / 4;
4097         int      index = (&hex->_i)[iDir];
4098         double linkLen = _grid->_coords[iDir][index+1] - _grid->_coords[iDir][index];
4099         if ( u - _grid->_tol > linkLen )
4100           continue;
4101
4102         if ( u < _grid->_tol ||
4103              u > linkLen - _grid->_tol ) // intersection at grid node
4104         {
4105           int  i = ! ( u < _grid->_tol ); // [0,1]
4106           int iN = link._nodes[ i ] - hex->_hexNodes; // [0-7]
4107
4108           const F_IntersectPoint * & ip = _grid->_gridIntP[ hex->_origNodeInd +
4109                                                             _grid->_nodeShift[iN] ];
4110           if ( !ip )
4111           {
4112             ip = _grid->_extIntPool.getNew();
4113             ip->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
4114             //ip->_transition = Trans_INTERNAL;
4115           }
4116           else if ( ip->_faceIDs.back() != _grid->PseudoIntExtFaceID() )
4117           {
4118             ip->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
4119           }
4120           hex->_nbFaceIntNodes++;
4121           //isCut = true;
4122         }
4123         else
4124         {
4125           const gp_Pnt&      p = intersection.Point( 1 );
4126           F_IntersectPoint* ip = _grid->_extIntPool.getNew();
4127           ip->_node = meshDS->AddNode( p.X(), p.Y(), p.Z() );
4128           ip->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
4129           ip->_transition = Trans_INTERNAL;
4130           meshDS->SetNodeInVolume( ip->_node, solidID );
4131
4132           CellsAroundLink fourCells( _grid, iDir );
4133           fourCells.Init( hex->_i, hex->_j, hex->_k, iLink );
4134           int i,j,k, cellIndex;
4135           for ( int iC = 0; iC < 4; ++iC ) // loop on 4 cells sharing the link
4136           {
4137             if ( !fourCells.GetCell( iC, i,j,k, cellIndex, iLink ))
4138               continue;
4139             Hexahedron * h = hexes[ cellIndex ];
4140             if ( !h )
4141               h = hexes[ cellIndex ] = new Hexahedron( *this, i, j, k, cellIndex );
4142             h->_hexLinks[iLink]._fIntPoints.push_back( ip );
4143             h->_nbFaceIntNodes++;
4144             //isCut = true;
4145           }
4146         }
4147       }
4148
4149       // if ( isCut )
4150       //   for ( size_t i = 0; i < hex->_eIntPoints.size(); ++i )
4151       //   {
4152       //     if ( _grid->IsInternal( hex->_eIntPoints[i]->_shapeID ) &&
4153       //          ! hex->_eIntPoints[i]->IsOnFace( _grid->PseudoIntExtFaceID() ))
4154       //       hex->_eIntPoints[i]->_faceIDs.push_back( _grid->PseudoIntExtFaceID() );
4155       //   }
4156       continue;
4157
4158     } // loop on all hexes
4159     return;
4160   }
4161
4162   //================================================================================
4163   /*!
4164    * \brief Return intersection point on INTERNAL FACE most distant from given ones
4165    */
4166   gp_Pnt Hexahedron::mostDistantInternalPnt( int hexIndex, const gp_Pnt& p1, const gp_Pnt& p2 )
4167   {
4168     gp_Pnt resultPnt = p1;
4169
4170     double maxDist2 = 0;
4171     for ( int iLink = 0; iLink < 12; ++iLink ) // check links
4172     {
4173       _Link& link = _hexLinks[ iLink ];
4174       for ( size_t i = 0; i < link._fIntPoints.size(); ++i )
4175         if ( _grid->PseudoIntExtFaceID() != link._fIntPoints[i]->_faceIDs[0] &&
4176              _grid->IsInternal( link._fIntPoints[i]->_faceIDs[0] ) &&
4177              link._fIntPoints[i]->_node )
4178         {
4179           gp_Pnt p = SMESH_NodeXYZ( link._fIntPoints[i]->_node );
4180           double d = p1.SquareDistance( p );
4181           if ( d > maxDist2 )
4182           {
4183             resultPnt = p;
4184             maxDist2  = d;
4185           }
4186           else
4187           {
4188             d = p2.SquareDistance( p );
4189             if ( d > maxDist2 )
4190             {
4191               resultPnt = p;
4192               maxDist2  = d;
4193             }
4194           }
4195         }
4196     }
4197     setIJK( hexIndex );
4198     _origNodeInd = _grid->NodeIndex( _i,_j,_k );
4199
4200     for ( size_t iN = 0; iN < 8; ++iN ) // check corners
4201     {
4202       _hexNodes[iN]._node     = _grid->_nodes   [ _origNodeInd + _grid->_nodeShift[iN] ];
4203       _hexNodes[iN]._intPoint = _grid->_gridIntP[ _origNodeInd + _grid->_nodeShift[iN] ];
4204       if ( _hexNodes[iN]._intPoint )
4205         for ( size_t iF = 0; iF < _hexNodes[iN]._intPoint->_faceIDs.size(); ++iF )
4206         {
4207           if ( _grid->IsInternal( _hexNodes[iN]._intPoint->_faceIDs[iF]))
4208           {
4209             gp_Pnt p = SMESH_NodeXYZ( _hexNodes[iN]._node );
4210             double d = p1.SquareDistance( p );
4211             if ( d > maxDist2 )
4212             {
4213               resultPnt = p;
4214               maxDist2  = d;
4215             }
4216             else
4217             {
4218               d = p2.SquareDistance( p );
4219               if ( d > maxDist2 )
4220               {
4221                 resultPnt = p;
4222                 maxDist2  = d;
4223               }
4224             }
4225           }
4226         }
4227     }
4228     if ( maxDist2 < _grid->_tol * _grid->_tol )
4229       return p1;
4230
4231     return resultPnt;
4232   }
4233
4234   //================================================================================
4235   /*!
4236    * \brief Finds intersection of a curve with a plane
4237    *  \param [in] u1 - parameter of one curve point
4238    *  \param [in] proj1 - projection of the curve point to the plane normal
4239    *  \param [in] u2 - parameter of another curve point
4240    *  \param [in] proj2 - projection of the other curve point to the plane normal
4241    *  \param [in] proj - projection of a point where the curve intersects the plane
4242    *  \param [in] curve - the curve
4243    *  \param [in] axis - the plane normal
4244    *  \param [in] origin - the plane origin
4245    *  \return gp_Pnt - the found intersection point
4246    */
4247   gp_Pnt Hexahedron::findIntPoint( double u1, double proj1,
4248                                    double u2, double proj2,
4249                                    double proj,
4250                                    BRepAdaptor_Curve& curve,
4251                                    const gp_XYZ& axis,
4252                                    const gp_XYZ& origin)
4253   {
4254     double r = (( proj - proj1 ) / ( proj2 - proj1 ));
4255     double u = u1 * ( 1 - r ) + u2 * r;
4256     gp_Pnt p = curve.Value( u );
4257     double newProj =  axis * ( p.XYZ() - origin );
4258     if ( Abs( proj - newProj ) > _grid->_tol / 10. )
4259     {
4260       if ( r > 0.5 )
4261         return findIntPoint( u2, proj2, u, newProj, proj, curve, axis, origin );
4262       else
4263         return findIntPoint( u1, proj2, u, newProj, proj, curve, axis, origin );
4264     }
4265     return p;
4266   }
4267
4268   //================================================================================
4269   /*!
4270    * \brief Returns indices of a hexahedron sub-entities holding a point
4271    *  \param [in] ip - intersection point
4272    *  \param [out] facets - 0-3 facets holding a point
4273    *  \param [out] sub - index of a vertex or an edge holding a point
4274    *  \return int - number of facets holding a point
4275    */
4276   int Hexahedron::getEntity( const E_IntersectPoint* ip, int* facets, int& sub )
4277   {
4278     enum { X = 1, Y = 2, Z = 4 }; // == 001, 010, 100
4279     int nbFacets = 0;
4280     int vertex = 0, edgeMask = 0;
4281
4282     if ( Abs( _grid->_coords[0][ _i   ] - ip->_uvw[0] ) < _grid->_tol ) {
4283       facets[ nbFacets++ ] = SMESH_Block::ID_F0yz;
4284       edgeMask |= X;
4285     }
4286     else if ( Abs( _grid->_coords[0][ _i+1 ] - ip->_uvw[0] ) < _grid->_tol ) {
4287       facets[ nbFacets++ ] = SMESH_Block::ID_F1yz;
4288       vertex   |= X;
4289       edgeMask |= X;
4290     }
4291     if ( Abs( _grid->_coords[1][ _j   ] - ip->_uvw[1] ) < _grid->_tol ) {
4292       facets[ nbFacets++ ] = SMESH_Block::ID_Fx0z;
4293       edgeMask |= Y;
4294     }
4295     else if ( Abs( _grid->_coords[1][ _j+1 ] - ip->_uvw[1] ) < _grid->_tol ) {
4296       facets[ nbFacets++ ] = SMESH_Block::ID_Fx1z;
4297       vertex   |= Y;
4298       edgeMask |= Y;
4299     }
4300     if ( Abs( _grid->_coords[2][ _k   ] - ip->_uvw[2] ) < _grid->_tol ) {
4301       facets[ nbFacets++ ] = SMESH_Block::ID_Fxy0;
4302       edgeMask |= Z;
4303     }
4304     else if ( Abs( _grid->_coords[2][ _k+1 ] - ip->_uvw[2] ) < _grid->_tol ) {
4305       facets[ nbFacets++ ] = SMESH_Block::ID_Fxy1;
4306       vertex   |= Z;
4307       edgeMask |= Z;
4308     }
4309
4310     switch ( nbFacets )
4311     {
4312     case 0: sub = 0;         break;
4313     case 1: sub = facets[0]; break;
4314     case 2: {
4315       const int edge [3][8] = {
4316         { SMESH_Block::ID_E00z, SMESH_Block::ID_E10z,
4317           SMESH_Block::ID_E01z, SMESH_Block::ID_E11z },
4318         { SMESH_Block::ID_E0y0, SMESH_Block::ID_E1y0, 0, 0,
4319           SMESH_Block::ID_E0y1, SMESH_Block::ID_E1y1 },
4320         { SMESH_Block::ID_Ex00, 0, SMESH_Block::ID_Ex10, 0,
4321           SMESH_Block::ID_Ex01, 0, SMESH_Block::ID_Ex11 }
4322       };
4323       switch ( edgeMask ) {
4324       case X | Y: sub = edge[ 0 ][ vertex ]; break;
4325       case X | Z: sub = edge[ 1 ][ vertex ]; break;
4326       default:    sub = edge[ 2 ][ vertex ];
4327       }
4328       break;
4329     }
4330     //case 3:
4331     default:
4332       sub = vertex + SMESH_Block::ID_FirstV;
4333     }
4334
4335     return nbFacets;
4336   }
4337   //================================================================================
4338   /*!
4339    * \brief Adds intersection with an EDGE
4340    */
4341   bool Hexahedron::addIntersection( const E_IntersectPoint* ip,
4342                                     vector< Hexahedron* >&  hexes,
4343                                     int ijk[], int dIJK[] )
4344   {
4345     bool added = false;
4346
4347     size_t hexIndex[4] = {
4348       _grid->CellIndex( ijk[0], ijk[1], ijk[2] ),
4349       dIJK[0] ? _grid->CellIndex( ijk[0]+dIJK[0], ijk[1], ijk[2] ) : -1,
4350       dIJK[1] ? _grid->CellIndex( ijk[0], ijk[1]+dIJK[1], ijk[2] ) : -1,
4351       dIJK[2] ? _grid->CellIndex( ijk[0], ijk[1], ijk[2]+dIJK[2] ) : -1
4352     };
4353     for ( int i = 0; i < 4; ++i )
4354     {
4355       if ( hexIndex[i] < hexes.size() && hexes[ hexIndex[i] ] )
4356       {
4357         Hexahedron* h = hexes[ hexIndex[i] ];
4358         h->_eIntPoints.reserve(2);
4359         h->_eIntPoints.push_back( ip );
4360         added = true;
4361
4362         // check if ip is really inside the hex
4363         if (SALOME::VerbosityActivated() && h->isOutParam( ip->_uvw ))
4364           throw SALOME_Exception("ip outside a hex");
4365       }
4366     }
4367     return added;
4368   }
4369   //================================================================================
4370   /*!
4371    * \brief Check if a hexahedron facet lies on a FACE
4372    *        Also return true if the facet does not interfere with any FACE
4373    */
4374   bool Hexahedron::isQuadOnFace( const size_t iQuad )
4375   {
4376     _Face& quad = _hexQuads[ iQuad ] ;
4377
4378     int nbGridNodesInt = 0; // nb FACE intersections at grid nodes
4379     int nbNoGeomNodes  = 0;
4380     for ( int iE = 0; iE < 4; ++iE )
4381     {
4382       nbNoGeomNodes = ( !quad._links[ iE ].FirstNode()->_intPoint &&
4383                         quad._links[ iE ].NbResultLinks() == 1      );
4384       nbGridNodesInt +=
4385         ( quad._links[ iE ].FirstNode()->_intPoint &&
4386           quad._links[ iE ].NbResultLinks() == 1   &&
4387           quad._links[ iE ].ResultLink( 0 ).FirstNode() == quad._links[ iE ].FirstNode() &&
4388           quad._links[ iE ].ResultLink( 0 ).LastNode()  == quad._links[ iE ].LastNode()   );
4389     }
4390     if ( nbNoGeomNodes == 4 )
4391       return true;
4392
4393     if ( nbGridNodesInt == 4 ) // all quad nodes are at FACE intersection
4394     {
4395       size_t iEmin = 0, minNbFaces = 1000;
4396       for ( int iE = 0; iE < 4; ++iE ) // look for a node with min nb FACEs
4397       {
4398         size_t nbFaces = quad._links[ iE ].FirstNode()->faces().size();
4399         if ( minNbFaces > nbFaces )
4400         {
4401           iEmin = iE;
4402           minNbFaces = nbFaces;
4403         }
4404       }
4405       // check if there is a FACE passing through all 4 nodes
4406       for ( const TGeomID& faceID : quad._links[ iEmin ].FirstNode()->faces() )
4407       {
4408         bool allNodesAtFace = true;
4409         for ( size_t iE = 0; iE < 4 &&  allNodesAtFace; ++iE )
4410           allNodesAtFace = ( iE == iEmin ||
4411                              quad._links[ iE ].FirstNode()->IsOnFace( faceID ));
4412         if ( allNodesAtFace ) // quad if on faceID
4413           return true;
4414       }
4415     }
4416     return false;
4417   }
4418   //================================================================================
4419   /*!
4420    * \brief Finds nodes at a path from one node to another via intersections with EDGEs
4421    */
4422   bool Hexahedron::findChain( _Node*          n1,
4423                               _Node*          n2,
4424                               _Face&          quad,
4425                               vector<_Node*>& chn )
4426   {
4427     chn.clear();
4428     chn.push_back( n1 );
4429     for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
4430       if ( !quad._eIntNodes[ iP ]->IsUsedInFace( &quad ) &&
4431            n1->IsLinked( quad._eIntNodes[ iP ]->_intPoint ) &&
4432            n2->IsLinked( quad._eIntNodes[ iP ]->_intPoint ))
4433       {
4434         chn.push_back( quad._eIntNodes[ iP ]);
4435         chn.push_back( n2 );
4436         quad._eIntNodes[ iP ]->_usedInFace = &quad;
4437         return true;
4438       }
4439     bool found;
4440     do
4441     {
4442       found = false;
4443       for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
4444         if ( !quad._eIntNodes[ iP ]->IsUsedInFace( &quad ) &&
4445              chn.back()->IsLinked( quad._eIntNodes[ iP ]->_intPoint ))
4446         {
4447           chn.push_back( quad._eIntNodes[ iP ]);
4448           found = ( quad._eIntNodes[ iP ]->_usedInFace = &quad );
4449           break;
4450         }
4451     } while ( found && ! chn.back()->IsLinked( n2->_intPoint ) );
4452
4453     if ( chn.back() != n2 && chn.back()->IsLinked( n2->_intPoint ))
4454       chn.push_back( n2 );
4455
4456     return chn.size() > 1;
4457   }
4458   //================================================================================
4459   /*!
4460    * \brief Try to heal a polygon whose ends are not connected
4461    */
4462   bool Hexahedron::closePolygon( _Face* polygon, vector<_Node*>& chainNodes ) const
4463   {
4464     int i = -1, nbLinks = polygon->_links.size();
4465     if ( nbLinks < 3 )
4466       return false;
4467     vector< _OrientedLink > newLinks;
4468     // find a node lying on the same FACE as the last one
4469     _Node*       node = polygon->_links.back().LastNode();
4470     TGeomID avoidFace = node->IsLinked( polygon->_links.back().FirstNode()->_intPoint );
4471     for ( i = nbLinks - 2; i >= 0; --i )
4472       if ( node->IsLinked( polygon->_links[i].FirstNode()->_intPoint, avoidFace ))
4473         break;
4474     if ( i >= 0 )
4475     {
4476       for ( ; i < nbLinks; ++i )
4477         newLinks.push_back( polygon->_links[i] );
4478     }
4479     else
4480     {
4481       // find a node lying on the same FACE as the first one
4482       node      = polygon->_links[0].FirstNode();
4483       avoidFace = node->IsLinked( polygon->_links[0].LastNode()->_intPoint );
4484       for ( i = 1; i < nbLinks; ++i )
4485         if ( node->IsLinked( polygon->_links[i].LastNode()->_intPoint, avoidFace ))
4486           break;
4487       if ( i < nbLinks )
4488         for ( nbLinks = i + 1, i = 0; i < nbLinks; ++i )
4489           newLinks.push_back( polygon->_links[i] );
4490     }
4491     if ( newLinks.size() > 1 )
4492     {
4493       polygon->_links.swap( newLinks );
4494       chainNodes.clear();
4495       chainNodes.push_back( polygon->_links.back().LastNode() );
4496       chainNodes.push_back( polygon->_links[0].FirstNode() );
4497       return true;
4498     }
4499     return false;
4500   }
4501   //================================================================================
4502   /*!
4503    * \brief Finds nodes on the same EDGE as the first node of avoidSplit.
4504    *
4505    * This function is for
4506    * 1) a case where an EDGE lies on a quad which lies on a FACE
4507    *    so that a part of quad in ON and another part is IN
4508    * 2) INTERNAL FACE passes through the 1st node of avoidSplit
4509    */
4510   bool Hexahedron::findChainOnEdge( const vector< _OrientedLink >& splits,
4511                                     const _OrientedLink&           prevSplit,
4512                                     const _OrientedLink&           avoidSplit,
4513                                     const std::set< TGeomID > &    concaveFaces,
4514                                     size_t &                       iS,
4515                                     _Face&                         quad,
4516                                     vector<_Node*>&                chn )
4517   {
4518     _Node* pn1 = prevSplit.FirstNode();
4519     _Node* pn2 = prevSplit.LastNode(); // pn2 is on EDGE, if not on INTERNAL FACE
4520     _Node* an3 = avoidSplit.LastNode();
4521     TGeomID avoidFace = pn1->IsLinked( pn2->_intPoint ); // FACE under the quad
4522     if ( avoidFace < 1 && pn1->_intPoint )
4523       return false;
4524
4525     chn.clear();
4526
4527     if ( !quad._eIntNodes.empty() ) // connect pn2 with EDGE intersections
4528     {
4529       chn.push_back( pn2 );
4530       bool found;
4531       do
4532       {
4533         found = false;
4534         for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
4535           if (( !quad._eIntNodes[ iP ]->IsUsedInFace( &quad )) &&
4536               ( chn.back()->IsLinked( quad._eIntNodes[ iP ]->_intPoint, avoidFace )) &&
4537               ( !avoidFace || quad._eIntNodes[ iP ]->IsOnFace( avoidFace )))
4538           {
4539             chn.push_back( quad._eIntNodes[ iP ]);
4540             found = ( quad._eIntNodes[ iP ]->_usedInFace = &quad );
4541             break;
4542           }
4543       } while ( found );
4544       pn2 = chn.back();
4545     }
4546
4547     _Node* n = 0, *stopNode = avoidSplit.LastNode();
4548
4549     if ( pn2 == prevSplit.LastNode() && // pn2 is at avoidSplit.FirstNode()
4550          !isCorner( stopNode ))         // stopNode is in the middle of a _hexLinks
4551     {
4552       // move stopNode to a _hexNodes
4553       for ( int iE = 0; iE < 4; ++iE ) // loop on 4 sides of a quadrangle
4554         for ( size_t iL = 0; iL < quad._links[ iE ].NbResultLinks(); ++iL )
4555         {
4556           const _Link* sideSplit = & quad._links[ iE ]._link->_splits[ iL ];
4557           if ( sideSplit == avoidSplit._link )
4558           {
4559             if ( quad._links[ iE ].LastNode()->Node() )
4560               stopNode = quad._links[ iE ].LastNode();
4561             iE = 4;
4562             break;
4563           }
4564         }
4565     }
4566
4567     // connect pn2 (probably new, at _eIntNodes) with a split
4568
4569     int i, iConn = 0;
4570     size_t nbCommon;
4571     TGeomID commonFaces[20];
4572     _Node* nPrev = nullptr;
4573     for ( i = splits.size()-1; i >= 0; --i )
4574     {
4575       if ( !splits[i] )
4576         continue;
4577
4578       bool stop = false;
4579       for ( int is1st = 0; is1st < 2; ++is1st )
4580       {
4581         _Node* nConn = is1st ? splits[i].FirstNode() : splits[i].LastNode();
4582         if ( nConn == nPrev )
4583         {
4584           if ( n == nConn )
4585             iConn = i;
4586           continue;
4587         }
4588         nPrev = nConn;
4589         if (( stop = ( nConn == stopNode )))
4590           break;
4591         // find a FACE connecting nConn with pn2 but not with an3
4592         if (( nConn != pn1 ) &&
4593             ( nConn->_intPoint && !nConn->_intPoint->_faceIDs.empty() ) &&
4594             ( nbCommon = nConn->GetCommonFaces( pn2->_intPoint, commonFaces )))
4595         {
4596           bool a3Coonect = true;
4597           for ( size_t iF = 0; iF < nbCommon &&  a3Coonect; ++iF )
4598             a3Coonect = an3->IsOnFace( commonFaces[ iF ]) || concaveFaces.count( commonFaces[ iF ]);
4599           if ( a3Coonect )
4600             continue;
4601
4602           if ( !n )
4603           {
4604             n = nConn;
4605             iConn = i + !is1st;
4606           }
4607           if ( nbCommon > 1 ) // nConn is linked with pn2 by an EDGE
4608           {
4609             n = nConn;
4610             iConn = i + !is1st;
4611             stop = true;
4612             break;
4613           }
4614         }
4615       }
4616       if ( stop )
4617       {
4618         i = iConn;
4619         break;
4620       }
4621     }
4622
4623     if ( n && n != stopNode )
4624     {
4625       if ( chn.empty() )
4626         chn.push_back( pn2 );
4627       chn.push_back( n );
4628       iS = i-1;
4629       return true;
4630     }
4631     else if ( !chn.empty() && chn.back()->_isInternalFlags )
4632     {
4633       // INTERNAL FACE partially cuts the quad
4634       for ( int ip = chn.size() - 2; ip >= 0; --ip )
4635         chn.push_back( chn[ ip ]);
4636       return true;
4637     }
4638     return false;
4639   }
4640   //================================================================================
4641   /*!
4642    * \brief Checks transition at the ginen intersection node of a link
4643    */
4644   bool Hexahedron::isOutPoint( _Link& link, int iP,
4645                                SMESH_MesherHelper& helper, const Solid* solid ) const
4646   {
4647     bool isOut = false;
4648
4649     if ( link._fIntNodes[iP]->faces().size() == 1 &&
4650          _grid->IsInternal( link._fIntNodes[iP]->face(0) ))
4651       return false;
4652
4653     const bool moreIntPoints = ( iP+1 < (int) link._fIntNodes.size() );
4654
4655     // get 2 _Node's
4656     _Node* n1 = link._fIntNodes[ iP ];
4657     if ( !n1->Node() )
4658       n1 = link._nodes[0];
4659     _Node* n2 = moreIntPoints ? link._fIntNodes[ iP+1 ] : 0;
4660     if ( !n2 || !n2->Node() )
4661       n2 = link._nodes[1];
4662     if ( !n2->Node() )
4663       return true;
4664
4665     // get all FACEs under n1 and n2
4666     set< TGeomID > faceIDs;
4667     if ( moreIntPoints ) faceIDs.insert( link._fIntNodes[iP+1]->faces().begin(),
4668                                          link._fIntNodes[iP+1]->faces().end() );
4669     if ( n2->_intPoint ) faceIDs.insert( n2->_intPoint->_faceIDs.begin(),
4670                                          n2->_intPoint->_faceIDs.end() );
4671     if ( faceIDs.empty() )
4672       return false; // n2 is inside
4673     if ( n1->_intPoint ) faceIDs.insert( n1->_intPoint->_faceIDs.begin(),
4674                                          n1->_intPoint->_faceIDs.end() );
4675     faceIDs.insert( link._fIntNodes[iP]->faces().begin(),
4676                     link._fIntNodes[iP]->faces().end() );
4677
4678     // get a point between 2 nodes
4679     gp_Pnt p1      = n1->Point();
4680     gp_Pnt p2      = n2->Point();
4681     gp_Pnt pOnLink = 0.8 * p1.XYZ() + 0.2 * p2.XYZ();
4682
4683     TopLoc_Location loc;
4684
4685     set< TGeomID >::iterator faceID = faceIDs.begin();
4686     for ( ; faceID != faceIDs.end(); ++faceID )
4687     {
4688       // project pOnLink on a FACE
4689       if ( *faceID < 1 || !solid->Contains( *faceID )) continue;
4690       const TopoDS_Face& face = TopoDS::Face( _grid->Shape( *faceID ));
4691       GeomAPI_ProjectPointOnSurf& proj = helper.GetProjector( face, loc, 0.1*_grid->_tol );
4692       gp_Pnt testPnt = pOnLink.Transformed( loc.Transformation().Inverted() );
4693       proj.Perform( testPnt );
4694       if ( proj.IsDone() && proj.NbPoints() > 0 )       
4695       {
4696         Standard_Real u,v;
4697         proj.LowerDistanceParameters( u,v );
4698
4699         if ( proj.LowerDistance() <= 0.1 * _grid->_tol )
4700         {
4701           isOut = false;
4702         }
4703         else
4704         {
4705           // find isOut by normals
4706           gp_Dir normal;
4707           if ( GeomLib::NormEstim( BRep_Tool::Surface( face, loc ),
4708                                    gp_Pnt2d( u,v ),
4709                                    0.1*_grid->_tol,
4710                                    normal ) < 3 )
4711           {
4712             if ( solid->Orientation( face ) == TopAbs_REVERSED )
4713               normal.Reverse();
4714             gp_Vec v( proj.NearestPoint(), testPnt );
4715             isOut = ( v * normal > 0 );
4716           }
4717         }
4718         if ( !isOut )
4719         {
4720           // classify a projection
4721           if ( !n1->IsOnFace( *faceID ) || !n2->IsOnFace( *faceID ))
4722           {
4723             BRepTopAdaptor_FClass2d cls( face, Precision::Confusion() );
4724             TopAbs_State state = cls.Perform( gp_Pnt2d( u,v ));
4725             if ( state == TopAbs_OUT )
4726             {
4727               isOut = true;
4728               continue;
4729             }
4730           }
4731           return false;
4732         }
4733       }
4734     }
4735     return isOut;
4736   }
4737   //================================================================================
4738   /*!
4739    * \brief Sort nodes on a FACE
4740    */
4741   void Hexahedron::sortVertexNodes(vector<_Node*>& nodes, _Node* curNode, TGeomID faceID)
4742   {
4743     if ( nodes.size() > 20 ) return;
4744
4745     // get shapes under nodes
4746     TGeomID nShapeIds[20], *nShapeIdsEnd = &nShapeIds[0] + nodes.size();
4747     for ( size_t i = 0; i < nodes.size(); ++i )
4748       if ( !( nShapeIds[i] = nodes[i]->ShapeID() ))
4749         return;
4750
4751     // get shapes of the FACE
4752     const TopoDS_Face&  face = TopoDS::Face( _grid->Shape( faceID ));
4753     list< TopoDS_Edge > edges;
4754     list< int >         nbEdges;
4755     int nbW = SMESH_Block::GetOrderedEdges (face, edges, nbEdges);
4756     if ( nbW > 1 ) {
4757       // select a WIRE - remove EDGEs of irrelevant WIREs from edges
4758       list< TopoDS_Edge >::iterator e = edges.begin(), eEnd = e;
4759       list< int >::iterator nE = nbEdges.begin();
4760       for ( ; nbW > 0; ++nE, --nbW )
4761       {
4762         std::advance( eEnd, *nE );
4763         for ( ; e != eEnd; ++e )
4764           for ( int i = 0; i < 2; ++i )
4765           {
4766             TGeomID id = i==0 ?
4767               _grid->ShapeID( *e ) :
4768               _grid->ShapeID( SMESH_MesherHelper::IthVertex( 0, *e ));
4769             if (( id > 0 ) &&
4770                 ( std::find( &nShapeIds[0], nShapeIdsEnd, id ) != nShapeIdsEnd ))
4771             {
4772               edges.erase( eEnd, edges.end() ); // remove rest wires
4773               e = eEnd = edges.end();
4774               --e;
4775               nbW = 0;
4776               break;
4777             }
4778           }
4779         if ( nbW > 0 )
4780           edges.erase( edges.begin(), eEnd ); // remove a current irrelevant wire
4781       }
4782     }
4783     // rotate edges to have the first one at least partially out of the hexa
4784     list< TopoDS_Edge >::iterator e = edges.begin(), eMidOut = edges.end();
4785     for ( ; e != edges.end(); ++e )
4786     {
4787       if ( !_grid->ShapeID( *e ))
4788         continue;
4789       bool isOut = false;
4790       gp_Pnt p;
4791       double uvw[3], f,l;
4792       for ( int i = 0; i < 2 && !isOut; ++i )
4793       {
4794         if ( i == 0 )
4795         {
4796           TopoDS_Vertex v = SMESH_MesherHelper::IthVertex( 0, *e );
4797           p = BRep_Tool::Pnt( v );
4798         }
4799         else if ( eMidOut == edges.end() )
4800         {
4801           TopLoc_Location loc;
4802           Handle(Geom_Curve) c = BRep_Tool::Curve( *e, loc, f, l);
4803           if ( c.IsNull() ) break;
4804           p = c->Value( 0.5 * ( f + l )).Transformed( loc );
4805         }
4806         else
4807         {
4808           continue;
4809         }
4810
4811         _grid->ComputeUVW( p.XYZ(), uvw );
4812         if ( isOutParam( uvw ))
4813         {
4814           if ( i == 0 )
4815             isOut = true;
4816           else
4817             eMidOut = e;
4818         }
4819       }
4820       if ( isOut )
4821         break;
4822     }
4823     if ( e != edges.end() )
4824       edges.splice( edges.end(), edges, edges.begin(), e );
4825     else if ( eMidOut != edges.end() )
4826       edges.splice( edges.end(), edges, edges.begin(), eMidOut );
4827
4828     // sort nodes according to the order of edges
4829     _Node*  orderNodes   [20];
4830     //TGeomID orderShapeIDs[20];
4831     size_t nbN = 0;
4832     TGeomID id, *pID = 0;
4833     for ( e = edges.begin(); e != edges.end(); ++e )
4834     {
4835       if (( id = _grid->ShapeID( SMESH_MesherHelper::IthVertex( 0, *e ))) &&
4836           (( pID = std::find( &nShapeIds[0], nShapeIdsEnd, id )) != nShapeIdsEnd ))
4837       {
4838         //orderShapeIDs[ nbN ] = id;
4839         orderNodes   [ nbN++ ] = nodes[ pID - &nShapeIds[0] ];
4840         *pID = -1;
4841       }
4842       if (( id = _grid->ShapeID( *e )) &&
4843           (( pID = std::find( &nShapeIds[0], nShapeIdsEnd, id )) != nShapeIdsEnd ))
4844       {
4845         //orderShapeIDs[ nbN ] = id;
4846         orderNodes   [ nbN++ ] = nodes[ pID - &nShapeIds[0] ];
4847         *pID = -1;
4848       }
4849     }
4850     if ( nbN != nodes.size() )
4851       return;
4852
4853     bool reverse = ( orderNodes[0    ]->Point().SquareDistance( curNode->Point() ) >
4854                      orderNodes[nbN-1]->Point().SquareDistance( curNode->Point() ));
4855
4856     for ( size_t i = 0; i < nodes.size(); ++i )
4857       nodes[ i ] = orderNodes[ reverse ? nbN-1-i : i ];
4858   }
4859
4860   //================================================================================
4861   /*!
4862    * \brief Adds computed elements to the mesh
4863    */
4864   int Hexahedron::addVolumes( SMESH_MesherHelper& helper )
4865   {
4866     F_IntersectPoint noIntPnt;
4867     const bool toCheckNodePos = _grid->IsToCheckNodePos();
4868     const bool useQuanta      = _grid->_toUseQuanta;
4869
4870     int nbAdded = 0;
4871     // add elements resulted from hexahedron intersection
4872     for ( _volumeDef* volDef = &_volumeDefs; volDef; volDef = volDef->_next )
4873     {
4874       vector< const SMDS_MeshNode* > nodes( volDef->_nodes.size() );
4875       for ( size_t iN = 0; iN < nodes.size(); ++iN )
4876       {
4877         if ( !( nodes[iN] = volDef->_nodes[iN].Node() ))
4878         {
4879           if ( const E_IntersectPoint* eip = volDef->_nodes[iN].EdgeIntPnt() )
4880           {
4881             nodes[iN] = volDef->_nodes[iN]._intPoint->_node =
4882               helper.AddNode( eip->_point.X(),
4883                               eip->_point.Y(),
4884                               eip->_point.Z() );
4885             if ( _grid->ShapeType( eip->_shapeID ) == TopAbs_VERTEX )
4886               helper.GetMeshDS()->SetNodeOnVertex( nodes[iN], eip->_shapeID );
4887             else
4888               helper.GetMeshDS()->SetNodeOnEdge( nodes[iN], eip->_shapeID );
4889           }
4890           else
4891             throw SALOME_Exception("Bug: no node at intersection point");
4892         }
4893         else if ( volDef->_nodes[iN]._intPoint &&
4894                   volDef->_nodes[iN]._intPoint->_node == volDef->_nodes[iN]._node )
4895         {
4896           // Update position of node at EDGE intersection;
4897           // see comment to _Node::Add( E_IntersectPoint )
4898           SMESHDS_Mesh* mesh = helper.GetMeshDS();
4899           TGeomID    shapeID = volDef->_nodes[iN].EdgeIntPnt()->_shapeID;
4900           mesh->UnSetNodeOnShape( nodes[iN] );
4901           if ( _grid->ShapeType( shapeID ) == TopAbs_VERTEX )
4902             mesh->SetNodeOnVertex( nodes[iN], shapeID );
4903           else
4904             mesh->SetNodeOnEdge( nodes[iN], shapeID );
4905         }
4906         else if ( toCheckNodePos &&
4907                   !nodes[iN]->isMarked() &&
4908                   _grid->ShapeType( nodes[iN]->GetShapeID() ) == TopAbs_FACE )
4909         {
4910           _grid->SetOnShape( nodes[iN], noIntPnt, /*v=*/nullptr,/*unset=*/true );
4911           nodes[iN]->setIsMarked( true );
4912         }
4913       } // loop to get nodes
4914
4915       const SMDS_MeshElement* v = 0;      
4916       if ( !volDef->_quantities.empty() )
4917       {                      
4918         if ( !useQuanta )
4919         {
4920           // split polyhedrons of with disjoint volumes
4921           std::vector<std::vector<int>> splitQuantities;
4922           std::vector<std::vector< const SMDS_MeshNode* > > splitNodes;
4923           if ( checkPolyhedronValidity( volDef, splitQuantities, splitNodes ) == 1 )
4924             v = addPolyhedronToMesh( volDef, helper, nodes, volDef->_quantities );
4925           else
4926           {
4927             int counter = -1;
4928             for (size_t id = 0; id < splitQuantities.size(); id++)
4929             {
4930               v = addPolyhedronToMesh( volDef, helper, splitNodes[ id ], splitQuantities[ id ] );
4931               if ( id < splitQuantities.size()-1 )
4932                 volDef->_brotherVolume.push_back( v );
4933               counter++;
4934             }
4935             nbAdded += counter;
4936           }
4937         }
4938         else
4939         {
4940           const double quanta = _grid->_quanta;
4941           double polyVol      = volDef->_size;
4942           double hexaVolume   = _sideLength[0] * _sideLength[1] * _sideLength[2];          
4943           if ( hexaVolume > 0.0 && polyVol/hexaVolume >= quanta /*set the volume if the relation is satisfied*/)
4944             v = helper.AddVolume( _hexNodes[0].BoundaryNode(), _hexNodes[2].BoundaryNode(),
4945                                   _hexNodes[3].BoundaryNode(), _hexNodes[1].BoundaryNode(),
4946                                   _hexNodes[4].BoundaryNode(), _hexNodes[6].BoundaryNode(),
4947                                   _hexNodes[7].BoundaryNode(), _hexNodes[5].BoundaryNode() );
4948           
4949         }
4950       }
4951       else
4952       {
4953         switch ( nodes.size() )
4954         {
4955         case 8: v = helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],
4956                                       nodes[4],nodes[5],nodes[6],nodes[7] );
4957           break;
4958         case 4: v = helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3] );
4959           break;
4960         case 6: v = helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],nodes[4],nodes[5] );
4961           break;
4962         case 5: v = helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],nodes[4] );
4963           break;
4964         }
4965       }
4966       volDef->_volume = v;
4967       nbAdded += bool( v );
4968
4969     } // loop on _volumeDefs chain
4970
4971     // avoid creating overlapping volumes (bos #24052)
4972     if ( nbAdded > 1 )
4973     {
4974       double sumSize = 0, maxSize = 0;
4975       _volumeDef* maxSizeDef = nullptr;
4976       for ( _volumeDef* volDef = &_volumeDefs; volDef; volDef = volDef->_next )
4977       {
4978         if ( !volDef->_volume )
4979           continue;
4980         sumSize += volDef->_size;
4981         if ( volDef->_size > maxSize )
4982         {
4983           maxSize    = volDef->_size;
4984           maxSizeDef = volDef;
4985         }
4986       }
4987       if ( sumSize > _sideLength[0] * _sideLength[1] * _sideLength[2] * 1.05 )
4988       {
4989         for ( _volumeDef* volDef = &_volumeDefs; volDef; volDef = volDef->_next )
4990           if ( volDef != maxSizeDef && volDef->_volume )
4991           {
4992             helper.GetMeshDS()->RemoveFreeElement( volDef->_volume, /*sm=*/nullptr,
4993                                                    /*fromGroups=*/false );
4994             volDef->_volume = nullptr;
4995             //volDef->_nodes.clear();
4996             --nbAdded;
4997           }
4998       }
4999     }
5000
5001     for ( _volumeDef* volDef = &_volumeDefs; volDef; volDef = volDef->_next )
5002     {
5003       if ( volDef->_volume )
5004       {
5005         helper.GetMeshDS()->SetMeshElementOnShape( volDef->_volume, volDef->_solidID );
5006         for (auto broVol : volDef->_brotherVolume )
5007         {
5008           helper.GetMeshDS()->SetMeshElementOnShape( broVol, volDef->_solidID );
5009         }
5010       }
5011     }
5012
5013     return nbAdded;
5014   }
5015   //================================================================================
5016   /*!
5017    * \brief Return true if the element is in a hole
5018    * \remark consider a cell to be in a hole if all links in any direction
5019    *          comes OUT of geometry
5020    */
5021   bool Hexahedron::isInHole() const
5022   {
5023     if ( !_vIntNodes.empty() )
5024       return false;
5025
5026     const size_t ijk[3] = { _i, _j, _k };
5027     F_IntersectPoint curIntPnt;
5028
5029     // consider a cell to be in a hole if all links in any direction
5030     // comes OUT of geometry
5031     for ( int iDir = 0; iDir < 3; ++iDir )
5032     {
5033       const vector<double>& coords = _grid->_coords[ iDir ];
5034       LineIndexer               li = _grid->GetLineIndexer( iDir );
5035       li.SetIJK( _i,_j,_k );
5036       size_t lineIndex[4] = { li.LineIndex  (),
5037                               li.LineIndex10(),
5038                               li.LineIndex01(),
5039                               li.LineIndex11() };
5040       bool allLinksOut = true, hasLinks = false;
5041       for ( int iL = 0; iL < 4 && allLinksOut; ++iL ) // loop on 4 links parallel to iDir
5042       {
5043         const _Link& link = _hexLinks[ iL + 4*iDir ];
5044         // check transition of the first node of a link
5045         const F_IntersectPoint* firstIntPnt = 0;
5046         if ( link._nodes[0]->Node() ) // 1st node is a hexa corner
5047         {
5048           curIntPnt._paramOnLine = coords[ ijk[ iDir ]] - coords[0] + _grid->_tol;
5049           const GridLine& line = _grid->_lines[ iDir ][ lineIndex[ iL ]];
5050           if ( !line._intPoints.empty() )
5051           {
5052             multiset< F_IntersectPoint >::const_iterator ip =
5053               line._intPoints.upper_bound( curIntPnt );
5054             --ip;
5055             firstIntPnt = &(*ip);
5056           }
5057         }
5058         else if ( !link._fIntPoints.empty() )
5059         {
5060           firstIntPnt = link._fIntPoints[0];
5061         }
5062
5063         if ( firstIntPnt )
5064         {
5065           hasLinks = true;
5066           allLinksOut = ( firstIntPnt->_transition == Trans_OUT &&
5067                           !_grid->IsShared( firstIntPnt->_faceIDs[0] ));
5068         }
5069       }
5070       if ( hasLinks && allLinksOut )
5071         return true;
5072     }
5073     return false;
5074   }
5075
5076   //================================================================================
5077   /*!
5078    * \brief Check if a polyherdon has an edge lying on EDGE shared by strange FACE
5079    *        that will be meshed by other algo
5080    */
5081   bool Hexahedron::hasStrangeEdge() const
5082   {
5083     if ( _eIntPoints.size() < 2 )
5084       return false;
5085
5086     TopTools_MapOfShape edges;
5087     for ( size_t i = 0; i < _eIntPoints.size(); ++i )
5088     {
5089       if ( !_grid->IsStrangeEdge( _eIntPoints[i]->_shapeID ))
5090         continue;
5091       const TopoDS_Shape& s = _grid->Shape( _eIntPoints[i]->_shapeID );
5092       if ( s.ShapeType() == TopAbs_EDGE )
5093       {
5094         if ( ! edges.Add( s ))
5095           return true; // an EDGE encounters twice
5096       }
5097       else
5098       {
5099         PShapeIteratorPtr edgeIt = _grid->_helper->GetAncestors( s,
5100                                                                  *_grid->_helper->GetMesh(),
5101                                                                  TopAbs_EDGE );
5102         while ( const TopoDS_Shape* edge = edgeIt->next() )
5103           if ( ! edges.Add( *edge ))
5104             return true; // an EDGE encounters twice
5105       }
5106     }
5107     return false;
5108   }
5109
5110   //================================================================================
5111   /*!
5112    * \brief Return true if a polyhedron passes _sizeThreshold criterion
5113    */
5114   bool Hexahedron::checkPolyhedronSize( bool cutByInternalFace, double & volume) const
5115   {
5116     volume = 0;
5117
5118     if ( cutByInternalFace && !_grid->_toUseThresholdForInternalFaces )
5119     {
5120       // check if any polygon fully lies on shared/internal FACEs
5121       for ( size_t iP = 0; iP < _polygons.size(); ++iP )
5122       {
5123         const _Face& polygon = _polygons[iP];
5124         if ( polygon._links.empty() )
5125           continue;
5126         bool allNodesInternal = true;
5127         for ( size_t iL = 0; iL < polygon._links.size() &&  allNodesInternal; ++iL )
5128         {
5129           _Node* n = polygon._links[ iL ].FirstNode();
5130           allNodesInternal = (( n->IsCutByInternal() ) ||
5131                               ( n->_intPoint && _grid->IsAnyShared( n->_intPoint->_faceIDs )));
5132         }
5133         if ( allNodesInternal )
5134           return true;
5135       }
5136     }
5137     for ( size_t iP = 0; iP < _polygons.size(); ++iP )
5138     {
5139       const _Face& polygon = _polygons[iP];
5140       if ( polygon._links.empty() )
5141         continue;
5142       gp_XYZ area (0,0,0);
5143       gp_XYZ p1 = polygon._links[ 0 ].FirstNode()->Point().XYZ();
5144       for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
5145       {
5146         gp_XYZ p2 = polygon._links[ iL ].LastNode()->Point().XYZ();
5147         area += p1 ^ p2;
5148         p1 = p2;
5149       }
5150       volume += p1 * area;
5151     }
5152     volume /= 6;
5153
5154     if ( this->hasStrangeEdge() && volume > 1e-13 )
5155       return true;
5156
5157     double initVolume = _sideLength[0] * _sideLength[1] * _sideLength[2];
5158
5159     return volume > initVolume / _grid->_sizeThreshold;
5160   }
5161
5162   //================================================================================
5163   /*!
5164    * \brief Check that all faces in polyhedron are connected so a unique volume is defined.
5165    *        We test that it is possible to go from any node to all nodes in the polyhedron.
5166    *        The set of nodes that can be visit within then defines a unique element.
5167    *        In case more than one polyhedron is detected. The function return the set of quantities and nodes defining separates elements.
5168    *        Reference to issue #bos[38521][EDF] Generate polyhedron with separate volume.
5169    */
5170   int Hexahedron::checkPolyhedronValidity( _volumeDef* volDef, std::vector<std::vector<int>>& splitQuantities, 
5171                                            std::vector<std::vector<const SMDS_MeshNode*>>& splitNodes )
5172   {  
5173     int mySet = 1;
5174     std::map<int,int> numberOfSets; // define set id with the number of faces associated!
5175     if ( !volDef->_quantities.empty() )
5176     {
5177       auto connectivity = volDef->_quantities;
5178       int accum = 0;
5179       std::vector<bool> allFaces( connectivity.size(), false );
5180       std::set<int> elementSet;
5181       allFaces[ 0 ] = true; // the first node below to the first face
5182       size_t connectedFaces = 1;
5183       // Start filling the set with the nodes of the first face
5184       splitQuantities.push_back( { connectivity[ 0 ] } );
5185       splitNodes.push_back( { volDef->_nodes[ 0 ].Node() } );
5186       elementSet.insert( volDef->_nodes[ 0 ].Node()->GetID() );
5187       for (int n = 1; n < connectivity[ 0 ]; n++)
5188       {
5189         elementSet.insert( volDef->_nodes[ n ].Node()->GetID() );
5190         splitNodes.back().push_back( volDef->_nodes[ n ].Node() );
5191       }
5192       
5193       numberOfSets.insert( std::pair<int,int>(mySet,1) );
5194       while ( connectedFaces != allFaces.size() )
5195       {
5196         for (size_t innerId = 1; innerId < connectivity.size(); innerId++)
5197         {
5198           if ( innerId == 1 )
5199             accum = connectivity[ 0 ];
5200
5201           if ( !allFaces[ innerId ] )
5202           {
5203             int faceCounter = 0;
5204             for (int n = 0; n < connectivity[ innerId ]; n++)
5205             {
5206               int nodeId = volDef->_nodes[ accum + n ].Node()->GetID();
5207               if ( elementSet.count( nodeId ) != 0 ) 
5208                 faceCounter++;
5209             }
5210             if ( faceCounter >= 2 ) // found coincidences nodes
5211             {
5212               for (int n = 0; n < connectivity[ innerId ]; n++)
5213               {
5214                 int nodeId = volDef->_nodes[ accum + n ].Node()->GetID();
5215                 // insert new nodes so other faces can be identified as belowing to the element
5216                 splitNodes.back().push_back( volDef->_nodes[ accum + n ].Node() );
5217                 elementSet.insert( nodeId );
5218               }
5219               allFaces[ innerId ] = true;
5220               splitQuantities.back().push_back( connectivity[ innerId ] );
5221               numberOfSets[ mySet ]++;
5222               connectedFaces++;
5223               innerId = 0; // to restart searching!
5224             }
5225           }
5226           accum += connectivity[ innerId ];
5227         }
5228
5229         if ( connectedFaces != allFaces.size() )
5230         {
5231           // empty the set, and fill it with nodes of a unvisited face!
5232           elementSet.clear();
5233           accum = connectivity[ 0 ];
5234           for (size_t faceId = 1; faceId < connectivity.size(); faceId++)
5235           {
5236             if ( !allFaces[ faceId ] )
5237             {
5238               splitNodes.push_back( { volDef->_nodes[ accum ].Node() } );
5239               elementSet.insert( volDef->_nodes[ accum ].Node()->GetID() );
5240               for (int n = 1; n < connectivity[ faceId ]; n++)
5241               {
5242                 elementSet.insert( volDef->_nodes[ accum + n ].Node()->GetID() );
5243                 splitNodes.back().push_back( volDef->_nodes[ accum + n ].Node() );
5244               }
5245
5246               splitQuantities.push_back( { connectivity[ faceId ] } );
5247               allFaces[ faceId ] = true;
5248               connectedFaces++;
5249               break;
5250             }
5251             accum += connectivity[ faceId ];
5252           }
5253           mySet++;
5254           numberOfSets.insert( std::pair<int,int>(mySet,1) );
5255         }
5256       }
5257
5258       if ( numberOfSets.size() > 1 )
5259       {
5260         bool allMoreThan2Faces = true;
5261         for( auto k : numberOfSets )
5262         {
5263           if ( k.second <= 2 )
5264             allMoreThan2Faces &= false;
5265         }
5266         
5267         if ( allMoreThan2Faces )
5268         {        
5269           // The separate objects are suspect to be closed
5270           return numberOfSets.size();        
5271         }
5272         else
5273         {
5274           // Have to index the last face nodes to the final set
5275           // contrary case return as it were a valid polyhedron for backward compatibility
5276           return 1;  
5277         }
5278       }
5279     }
5280     return numberOfSets.size();
5281   }
5282
5283   
5284   //================================================================================
5285   /*!
5286    * \brief add original or separated polyhedrons to the mesh
5287    */
5288   const SMDS_MeshElement* Hexahedron::addPolyhedronToMesh( _volumeDef* volDef,  SMESH_MesherHelper& helper, const std::vector<const SMDS_MeshNode*>& nodes, 
5289                                                             const std::vector<int>& quantities )
5290   {
5291     const SMDS_MeshElement* v = helper.AddPolyhedralVolume( nodes, quantities );
5292
5293     volDef->_size = SMDS_VolumeTool( v ).GetSize();
5294     if ( volDef->_size < 0 ) // invalid polyhedron
5295     {
5296       if ( ! SMESH_MeshEditor( helper.GetMesh() ).Reorient( v ) || // try to fix
5297           SMDS_VolumeTool( v ).GetSize() < 0 )
5298       {
5299         helper.GetMeshDS()->RemoveFreeElement( v, /*sm=*/nullptr, /*fromGroups=*/false );
5300         v = nullptr;
5301         //_hasTooSmall = true;
5302
5303         if (SALOME::VerbosityActivated())
5304         {
5305           std::cout << "Remove INVALID polyhedron, _cellID = " << _cellID
5306                     << " ijk = ( " << _i << " " << _j << " " << _k << " ) "
5307                     << " solid " << volDef->_solidID << std::endl;
5308         }
5309       }
5310     }
5311     return v;
5312   }
5313
5314   //================================================================================
5315   /*!
5316    * \brief Tries to create a hexahedron
5317    */
5318   bool Hexahedron::addHexa()
5319   {
5320     int nbQuad = 0, iQuad = -1;
5321     for ( size_t i = 0; i < _polygons.size(); ++i )
5322     {
5323       if ( _polygons[i]._links.empty() )
5324         continue;
5325       if ( _polygons[i]._links.size() != 4 )
5326         return false;
5327       ++nbQuad;
5328       if ( iQuad < 0 )
5329         iQuad = i;
5330     }
5331     if ( nbQuad != 6 )
5332       return false;
5333
5334     _Node* nodes[8];
5335     int nbN = 0;
5336     for ( int iL = 0; iL < 4; ++iL )
5337     {
5338       // a base node
5339       nodes[iL] = _polygons[iQuad]._links[iL].FirstNode();
5340       ++nbN;
5341
5342       // find a top node above the base node
5343       _Link* link = _polygons[iQuad]._links[iL]._link;
5344       if ( !link->_faces[0] || !link->_faces[1] )
5345         return debugDumpLink( link );
5346       // a quadrangle sharing <link> with _polygons[iQuad]
5347       _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[iQuad] )];
5348       for ( int i = 0; i < 4; ++i )
5349         if ( quad->_links[i]._link == link )
5350         {
5351           // 1st node of a link opposite to <link> in <quad>
5352           nodes[iL+4] = quad->_links[(i+2)%4].FirstNode();
5353           ++nbN;
5354           break;
5355         }
5356     }
5357     if ( nbN == 8 )
5358       _volumeDefs.Set( &nodes[0], 8 );
5359
5360     return nbN == 8;
5361   }
5362   //================================================================================
5363   /*!
5364    * \brief Tries to create a tetrahedron
5365    */
5366   bool Hexahedron::addTetra()
5367   {
5368     int iTria = -1;
5369     for ( size_t i = 0; i < _polygons.size() && iTria < 0; ++i )
5370       if ( _polygons[i]._links.size() == 3 )
5371         iTria = i;
5372     if ( iTria < 0 )
5373       return false;
5374
5375     _Node* nodes[4];
5376     nodes[0] = _polygons[iTria]._links[0].FirstNode();
5377     nodes[1] = _polygons[iTria]._links[1].FirstNode();
5378     nodes[2] = _polygons[iTria]._links[2].FirstNode();
5379
5380     _Link* link = _polygons[iTria]._links[0]._link;
5381     if ( !link->_faces[0] || !link->_faces[1] )
5382       return debugDumpLink( link );
5383
5384     // a triangle sharing <link> with _polygons[0]
5385     _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[iTria] )];
5386     for ( int i = 0; i < 3; ++i )
5387       if ( tria->_links[i]._link == link )
5388       {
5389         nodes[3] = tria->_links[(i+1)%3].LastNode();
5390         _volumeDefs.Set( &nodes[0], 4 );
5391         return true;
5392       }
5393
5394     return false;
5395   }
5396   //================================================================================
5397   /*!
5398    * \brief Tries to create a pentahedron
5399    */
5400   bool Hexahedron::addPenta()
5401   {
5402     // find a base triangular face
5403     int iTri = -1;
5404     for ( int iF = 0; iF < 5 && iTri < 0; ++iF )
5405       if ( _polygons[ iF ]._links.size() == 3 )
5406         iTri = iF;
5407     if ( iTri < 0 ) return false;
5408
5409     // find nodes
5410     _Node* nodes[6];
5411     int nbN = 0;
5412     for ( int iL = 0; iL < 3; ++iL )
5413     {
5414       // a base node
5415       nodes[iL] = _polygons[ iTri ]._links[iL].FirstNode();
5416       ++nbN;
5417
5418       // find a top node above the base node
5419       _Link* link = _polygons[ iTri ]._links[iL]._link;
5420       if ( !link->_faces[0] || !link->_faces[1] )
5421         return debugDumpLink( link );
5422       // a quadrangle sharing <link> with a base triangle
5423       _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[ iTri ] )];
5424       if ( quad->_links.size() != 4 ) return false;
5425       for ( int i = 0; i < 4; ++i )
5426         if ( quad->_links[i]._link == link )
5427         {
5428           // 1st node of a link opposite to <link> in <quad>
5429           nodes[iL+3] = quad->_links[(i+2)%4].FirstNode();
5430           ++nbN;
5431           break;
5432         }
5433     }
5434     if ( nbN == 6 )
5435       _volumeDefs.Set( &nodes[0], 6 );
5436
5437     return ( nbN == 6 );
5438   }
5439   //================================================================================
5440   /*!
5441    * \brief Tries to create a pyramid
5442    */
5443   bool Hexahedron::addPyra()
5444   {
5445     // find a base quadrangle
5446     int iQuad = -1;
5447     for ( int iF = 0; iF < 5 && iQuad < 0; ++iF )
5448       if ( _polygons[ iF ]._links.size() == 4 )
5449         iQuad = iF;
5450     if ( iQuad < 0 ) return false;
5451
5452     // find nodes
5453     _Node* nodes[5];
5454     nodes[0] = _polygons[iQuad]._links[0].FirstNode();
5455     nodes[1] = _polygons[iQuad]._links[1].FirstNode();
5456     nodes[2] = _polygons[iQuad]._links[2].FirstNode();
5457     nodes[3] = _polygons[iQuad]._links[3].FirstNode();
5458
5459     _Link* link = _polygons[iQuad]._links[0]._link;
5460     if ( !link->_faces[0] || !link->_faces[1] )
5461       return debugDumpLink( link );
5462
5463     // a triangle sharing <link> with a base quadrangle
5464     _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[ iQuad ] )];
5465     if ( tria->_links.size() != 3 ) return false;
5466     for ( int i = 0; i < 3; ++i )
5467       if ( tria->_links[i]._link == link )
5468       {
5469         nodes[4] = tria->_links[(i+1)%3].LastNode();
5470         _volumeDefs.Set( &nodes[0], 5 );
5471         return true;
5472       }
5473
5474     return false;
5475   }
5476   //================================================================================
5477   /*!
5478    * \brief Return true if there are _eIntPoints at EDGEs forming a  concave corner
5479    */
5480   bool Hexahedron::hasEdgesAround( const ConcaveFace* cf ) const
5481   {
5482     int nbEdges = 0;
5483     ConcaveFace foundGeomHolder;
5484     for ( const E_IntersectPoint* ip : _eIntPoints )
5485     {
5486       if ( cf->HasEdge( ip->_shapeID ))
5487       {
5488         if ( ++nbEdges == 2 )
5489           return true;
5490         foundGeomHolder.SetEdge( ip->_shapeID );
5491       }
5492       else if ( ip->_faceIDs.size() >= 3 )
5493       {
5494         const TGeomID & vID = ip->_shapeID;
5495         if ( cf->HasVertex( vID ) && !foundGeomHolder.HasVertex( vID ))
5496         {
5497           if ( ++nbEdges == 2 )
5498             return true;
5499           foundGeomHolder.SetVertex( vID );
5500         }
5501       }
5502     }
5503
5504     for ( const _Node& hexNode: _hexNodes )
5505     {
5506       if ( !hexNode._node || !hexNode._intPoint )
5507         continue;
5508       const B_IntersectPoint* ip = hexNode._intPoint;
5509       if ( ip->_faceIDs.size() == 2 ) // EDGE
5510       {
5511         TGeomID edgeID = hexNode._node->GetShapeID();
5512         if ( cf->HasEdge( edgeID ) && !foundGeomHolder.HasEdge( edgeID ))
5513         {
5514           foundGeomHolder.SetEdge( edgeID );
5515           if ( ++nbEdges == 2 )
5516             return true;
5517         }
5518       }
5519       else if ( ip->_faceIDs.size() >= 3 ) // VERTEX
5520       {
5521         TGeomID vID = hexNode._node->GetShapeID();
5522         if ( cf->HasVertex( vID ) && !foundGeomHolder.HasVertex( vID ))
5523         {
5524           if ( ++nbEdges == 2 )
5525             return true;
5526           foundGeomHolder.SetVertex( vID );
5527         }
5528       }
5529     }
5530
5531     return false;
5532   }
5533   //================================================================================
5534   /*!
5535    * \brief Dump a link and return \c false
5536    */
5537   bool Hexahedron::debugDumpLink( Hexahedron::_Link* link )
5538   {
5539     if (SALOME::VerbosityActivated())
5540     {
5541       gp_Pnt p1 = link->_nodes[0]->Point(), p2 = link->_nodes[1]->Point();
5542       cout << "BUG: not shared link. IKJ = ( "<< _i << " " << _j << " " << _k << " )" << endl
5543           << "n1 (" << p1.X() << ", "<< p1.Y() << ", "<< p1.Z() << " )" << endl
5544           << "n2 (" << p2.X() << ", "<< p2.Y() << ", "<< p2.Z() << " )" << endl;
5545     }
5546
5547     return false;
5548   }
5549   //================================================================================
5550   /*!
5551    * \brief Classify a point by grid parameters
5552    */
5553   bool Hexahedron::isOutParam(const double uvw[3]) const
5554   {
5555     return (( _grid->_coords[0][ _i   ] - _grid->_tol > uvw[0] ) ||
5556             ( _grid->_coords[0][ _i+1 ] + _grid->_tol < uvw[0] ) ||
5557             ( _grid->_coords[1][ _j   ] - _grid->_tol > uvw[1] ) ||
5558             ( _grid->_coords[1][ _j+1 ] + _grid->_tol < uvw[1] ) ||
5559             ( _grid->_coords[2][ _k   ] - _grid->_tol > uvw[2] ) ||
5560             ( _grid->_coords[2][ _k+1 ] + _grid->_tol < uvw[2] ));
5561   }
5562   //================================================================================
5563   /*!
5564    * \brief Find existing triangulation of a polygon
5565    */
5566   int findExistingTriangulation( const SMDS_MeshElement*              polygon,
5567                                  //const SMDS_Mesh*                     mesh,
5568                                  std::vector< const SMDS_MeshNode* >& nodes )
5569   {
5570     int nbSplits = 0;
5571     nodes.clear();
5572     std::vector<const SMDS_MeshNode *>    twoNodes(2);
5573     std::vector<const SMDS_MeshElement *> foundFaces; foundFaces.reserve(10);
5574     std::set< const SMDS_MeshElement * >  avoidFaces; avoidFaces.insert( polygon );
5575
5576     const int nbPolyNodes = polygon->NbCornerNodes();
5577     twoNodes[1] = polygon->GetNode( nbPolyNodes - 1 );
5578     for ( int iN = 0; iN < nbPolyNodes; ++iN ) // loop on border links of polygon
5579     {
5580       twoNodes[0] = polygon->GetNode( iN );
5581
5582       int nbFaces = SMDS_Mesh::GetElementsByNodes( twoNodes, foundFaces, SMDSAbs_Face );
5583       int nbOkFaces = 0;
5584       for ( int iF = 0; iF < nbFaces; ++iF ) // keep faces lying over polygon
5585       {
5586         if ( avoidFaces.count( foundFaces[ iF ]))
5587           continue;
5588         int i, nbFaceNodes = foundFaces[ iF ]->NbCornerNodes();
5589         for ( i = 0; i < nbFaceNodes; ++i )
5590         {
5591           const SMDS_MeshNode* n = foundFaces[ iF ]->GetNode( i );
5592           bool isCommonNode = ( n == twoNodes[0] ||
5593                                 n == twoNodes[1] ||
5594                                 polygon->GetNodeIndex( n ) >= 0 );
5595           if ( !isCommonNode )
5596             break;
5597         }
5598         if ( i == nbFaceNodes ) // all nodes of foundFaces[iF] are shared with polygon
5599           if ( nbOkFaces++ != iF )
5600             foundFaces[ nbOkFaces-1 ] = foundFaces[ iF ];
5601       }
5602       if ( nbOkFaces > 0 )
5603       {
5604         int iFaceSelected = 0;
5605         if ( nbOkFaces > 1 ) // select a face with minimal distance from polygon
5606         {
5607           double minDist = Precision::Infinite();
5608           for ( int iF = 0; iF < nbOkFaces; ++iF )
5609           {
5610             int i, nbFaceNodes = foundFaces[ iF ]->NbCornerNodes();
5611             gp_XYZ gc = SMESH_NodeXYZ( foundFaces[ iF ]->GetNode( 0 ));
5612             for ( i = 1; i < nbFaceNodes; ++i )
5613               gc += SMESH_NodeXYZ( foundFaces[ iF ]->GetNode( i ));
5614             gc /= nbFaceNodes;
5615
5616             double dist = SMESH_MeshAlgos::GetDistance( polygon, gc );
5617             if ( dist < minDist )
5618             {
5619               minDist = dist;
5620               iFaceSelected = iF;
5621             }
5622           }
5623         }
5624         if ( foundFaces[ iFaceSelected ]->NbCornerNodes() != 3 )
5625           return 0;
5626         nodes.insert( nodes.end(),
5627                       foundFaces[ iFaceSelected ]->begin_nodes(),
5628                       foundFaces[ iFaceSelected ]->end_nodes());
5629         if ( !SMESH_MeshAlgos::IsRightOrder( foundFaces[ iFaceSelected ],
5630                                              twoNodes[0], twoNodes[1] ))
5631         {
5632           // reverse just added nodes
5633           std::reverse( nodes.end() - 3, nodes.end() );
5634         }
5635         avoidFaces.insert( foundFaces[ iFaceSelected ]);
5636         nbSplits++;
5637       }
5638
5639       twoNodes[1] = twoNodes[0];
5640
5641     } // loop on polygon nodes
5642
5643     return nbSplits;
5644   }
5645   //================================================================================
5646   /*!
5647    * \brief Divide a polygon into triangles and modify accordingly an adjacent polyhedron
5648    */
5649   void splitPolygon( const SMDS_MeshElement*         polygon,
5650                      SMDS_VolumeTool &               volume,
5651                      const int                       facetIndex,
5652                      const TGeomID                   faceID,
5653                      const TGeomID                   solidID,
5654                      SMESH_MeshEditor::ElemFeatures& face,
5655                      SMESH_MeshEditor&               editor,
5656                      const bool                      reinitVolume)
5657   {
5658     SMESH_MeshAlgos::Triangulate divider(/*optimize=*/false);
5659     bool triangulationExist = false;
5660     int nbTrias = findExistingTriangulation( polygon, face.myNodes );
5661     if ( nbTrias > 0 )
5662       triangulationExist = true;
5663     else
5664       nbTrias = divider.GetTriangles( polygon, face.myNodes );
5665     face.myNodes.resize( nbTrias * 3 );
5666
5667     SMESH_MeshEditor::ElemFeatures newVolumeDef;
5668     newVolumeDef.Init( volume.Element() );
5669     newVolumeDef.SetID( volume.Element()->GetID() );
5670
5671     newVolumeDef.myPolyhedQuantities.reserve( volume.NbFaces() + nbTrias );
5672     newVolumeDef.myNodes.reserve( volume.NbNodes() + nbTrias * 3 );
5673
5674     SMESHDS_Mesh* meshDS = editor.GetMeshDS();
5675     SMDS_MeshElement* newTriangle;
5676     for ( int iF = 0, nF = volume.NbFaces(); iF < nF; iF++ )
5677     {
5678       if ( iF == facetIndex )
5679       {
5680         newVolumeDef.myPolyhedQuantities.push_back( 3 );
5681         newVolumeDef.myNodes.insert( newVolumeDef.myNodes.end(),
5682                                      face.myNodes.begin(),
5683                                      face.myNodes.begin() + 3 );
5684         meshDS->RemoveFreeElement( polygon, 0, false );
5685         if ( !triangulationExist )
5686         {
5687           newTriangle = meshDS->AddFace( face.myNodes[0], face.myNodes[1], face.myNodes[2] );
5688           meshDS->SetMeshElementOnShape( newTriangle, faceID );
5689         }
5690       }
5691       else
5692       {
5693         const SMDS_MeshNode** nn = volume.GetFaceNodes( iF );
5694         const size_t nbFaceNodes = volume.NbFaceNodes ( iF );
5695         newVolumeDef.myPolyhedQuantities.push_back( nbFaceNodes );
5696         newVolumeDef.myNodes.insert( newVolumeDef.myNodes.end(), nn, nn + nbFaceNodes );
5697       }
5698     }
5699
5700     for ( size_t iN = 3; iN < face.myNodes.size(); iN += 3 )
5701     {
5702       newVolumeDef.myPolyhedQuantities.push_back( 3 );
5703       newVolumeDef.myNodes.insert( newVolumeDef.myNodes.end(),
5704                                    face.myNodes.begin() + iN,
5705                                    face.myNodes.begin() + iN + 3 );
5706       if ( !triangulationExist )
5707       {
5708         newTriangle = meshDS->AddFace( face.myNodes[iN], face.myNodes[iN+1], face.myNodes[iN+2] );
5709         meshDS->SetMeshElementOnShape( newTriangle, faceID );
5710       }
5711     }
5712
5713     meshDS->RemoveFreeElement( volume.Element(), 0, false );
5714     SMDS_MeshElement* newVolume = editor.AddElement( newVolumeDef.myNodes, newVolumeDef );
5715     meshDS->SetMeshElementOnShape( newVolume, solidID );
5716
5717     if ( reinitVolume )
5718     {
5719       volume.Set( 0 );
5720       volume.Set( newVolume );
5721     }
5722     return;
5723   }
5724   //================================================================================
5725   /*!
5726    * \brief Look for a FACE supporting all given nodes made on EDGEs and VERTEXes
5727    */
5728   TGeomID findCommonFace( const std::vector< const SMDS_MeshNode* > & nn,
5729                           const SMESH_Mesh*                           mesh )
5730   {
5731     TGeomID faceID = 0;
5732     TGeomID shapeIDs[20];
5733     for ( size_t iN = 0; iN < nn.size(); ++iN )
5734       shapeIDs[ iN ] = nn[ iN ]->GetShapeID();
5735
5736     SMESH_subMesh* sm = mesh->GetSubMeshContaining( shapeIDs[ 0 ]);
5737     for ( const SMESH_subMesh * smFace : sm->GetAncestors() )
5738     {
5739       if ( smFace->GetSubShape().ShapeType() != TopAbs_FACE )
5740         continue;
5741
5742       faceID = smFace->GetId();
5743
5744       for ( size_t iN = 1; iN < nn.size() &&  faceID; ++iN )
5745       {
5746         if ( !smFace->DependsOn( shapeIDs[ iN ]))
5747           faceID = 0;
5748       }
5749       if ( faceID > 0 )
5750         break;
5751     }
5752     return faceID;
5753   }
5754   //================================================================================
5755   /*!
5756    * \brief Create mesh faces at free facets
5757    */
5758   void Hexahedron::addFaces( SMESH_MesherHelper&                       helper,
5759                              const vector< const SMDS_MeshElement* > & boundaryVolumes )
5760   {
5761     if ( !_grid->_toCreateFaces )
5762       return;
5763
5764     SMDS_VolumeTool vTool;
5765     vector<int> bndFacets;
5766     SMESH_MeshEditor editor( helper.GetMesh() );
5767     SMESH_MeshEditor::ElemFeatures face( SMDSAbs_Face );
5768     SMESHDS_Mesh* meshDS = helper.GetMeshDS();
5769
5770     bool isQuantaSet =  _grid->_toUseQuanta;    
5771     // check if there are internal or shared FACEs
5772     bool hasInternal = ( !_grid->_geometry.IsOneSolid() ||
5773                          _grid->_geometry._soleSolid.HasInternalFaces() );   
5774
5775     for ( size_t iV = 0; iV < boundaryVolumes.size(); ++iV )
5776     {
5777       if ( !vTool.Set( boundaryVolumes[ iV ]))
5778         continue;
5779       TGeomID solidID = vTool.Element()->GetShapeID();
5780       Solid *   solid = _grid->GetOneOfSolids( solidID );
5781
5782       // find boundary facets
5783       bndFacets.clear();
5784       for ( int iF = 0, n = vTool.NbFaces(); iF < n; iF++ )
5785       {
5786         const SMDS_MeshElement* otherVol;
5787         bool isBoundary = isQuantaSet ? vTool.IsFreeFaceCheckAllNodes( iF, &otherVol ) : vTool.IsFreeFace( iF, &otherVol );
5788         if ( isBoundary )
5789         {
5790           bndFacets.push_back( iF );
5791         }
5792         else if (( hasInternal ) ||
5793                  ( !_grid->IsSolid( otherVol->GetShapeID() )))
5794         {
5795           // check if all nodes are on internal/shared FACEs
5796           isBoundary = true;
5797           const SMDS_MeshNode** nn = vTool.GetFaceNodes( iF );
5798           const size_t nbFaceNodes = vTool.NbFaceNodes ( iF );
5799           for ( size_t iN = 0; iN < nbFaceNodes &&  isBoundary; ++iN )
5800             isBoundary = ( nn[ iN ]->GetShapeID() != solidID );
5801           if ( isBoundary )
5802             bndFacets.push_back( -( iF+1 )); // !!! minus ==> to check the FACE
5803         }
5804       }
5805       if ( bndFacets.empty() )
5806         continue;
5807
5808       // create faces
5809       if ( !vTool.IsPoly() )
5810         vTool.SetExternalNormal();
5811       for ( size_t i = 0; i < bndFacets.size(); ++i ) // loop on boundary facets
5812       {
5813         const bool    isBoundary = ( bndFacets[i] >= 0 );
5814         const int         iFacet = isBoundary ? bndFacets[i] : -bndFacets[i]-1;
5815         const SMDS_MeshNode** nn = vTool.GetFaceNodes( iFacet );
5816         const size_t nbFaceNodes = vTool.NbFaceNodes ( iFacet );
5817         face.myNodes.assign( nn, nn + nbFaceNodes );
5818
5819         TGeomID faceID = 0;
5820         const SMDS_MeshElement* existFace = 0, *newFace = 0;
5821
5822         if (( existFace = meshDS->FindElement( face.myNodes, SMDSAbs_Face )))
5823         {
5824           if ( existFace->isMarked() )
5825             continue; // created by this method
5826           faceID = existFace->GetShapeID();
5827         }
5828         else
5829         {
5830           // look for a supporting FACE
5831           for ( size_t iN = 0; iN < nbFaceNodes &&  !faceID; ++iN ) // look for a node on FACE
5832           {
5833             if ( nn[ iN ]->GetPosition()->GetDim() == 2 )
5834               faceID = nn[ iN ]->GetShapeID();
5835           }
5836           if ( faceID == 0 && !isQuantaSet /*if quanta is set boundary nodes at boundary does not coincide with any geometrical face */ )
5837             faceID = findCommonFace( face.myNodes, helper.GetMesh() );
5838
5839           bool toCheckFace = faceID && (( !isBoundary ) ||
5840                                         ( hasInternal && _grid->_toUseThresholdForInternalFaces ));
5841           if ( toCheckFace ) // check if all nodes are on the found FACE
5842           {
5843             SMESH_subMesh* faceSM = helper.GetMesh()->GetSubMeshContaining( faceID );
5844             for ( size_t iN = 0; iN < nbFaceNodes &&  faceID; ++iN )
5845             {
5846               TGeomID subID = nn[ iN ]->GetShapeID();
5847               if ( subID != faceID && !faceSM->DependsOn( subID ))
5848                 faceID = 0;
5849             }
5850             // if ( !faceID && !isBoundary )
5851             //   continue;
5852           }
5853           if ( !faceID && !isBoundary && !isQuantaSet )
5854             continue;
5855         }
5856
5857         // orient a new face according to supporting FACE orientation in shape_to_mesh
5858         if ( !isBoundary && !solid->IsOutsideOriented( faceID ))
5859         {
5860           if ( existFace )
5861             editor.Reorient( existFace );
5862           else
5863             std::reverse( face.myNodes.begin(), face.myNodes.end() );
5864         }
5865
5866         if ( ! ( newFace = existFace ))
5867         {
5868           face.SetPoly( nbFaceNodes > 4 );
5869           newFace = editor.AddElement( face.myNodes, face );
5870           if ( !newFace )
5871             continue;
5872           newFace->setIsMarked( true ); // to distinguish from face created in getBoundaryElems()
5873         }
5874
5875         if ( faceID && _grid->IsBoundaryFace( faceID )) // face is not shared
5876         {
5877           // set newFace to the found FACE provided that it fully lies on the FACE
5878           for ( size_t iN = 0; iN < nbFaceNodes &&  faceID; ++iN )
5879             if ( nn[iN]->GetShapeID() == solidID )
5880             {
5881               if ( existFace )
5882                 meshDS->UnSetMeshElementOnShape( existFace, _grid->Shape( faceID ));
5883               faceID = 0;
5884             }
5885         }
5886
5887         if ( faceID && nbFaceNodes > 4 &&
5888              !_grid->IsInternal( faceID ) &&
5889              !_grid->IsShared( faceID ) &&
5890              !_grid->IsBoundaryFace( faceID ))
5891         {
5892           // split a polygon that will be used by other 3D algorithm
5893           if ( !existFace )
5894             splitPolygon( newFace, vTool, iFacet, faceID, solidID,
5895                           face, editor, i+1 < bndFacets.size() );
5896         }
5897         else
5898         {
5899           if ( faceID )
5900             meshDS->SetMeshElementOnShape( newFace, faceID );
5901           else
5902             meshDS->SetMeshElementOnShape( newFace, solidID );
5903         }
5904       } // loop on bndFacets
5905     } // loop on boundaryVolumes
5906
5907
5908     // Orient coherently mesh faces on INTERNAL FACEs
5909
5910     if ( hasInternal )
5911     {
5912       TopExp_Explorer faceExp( _grid->_geometry._mainShape, TopAbs_FACE );
5913       for ( ; faceExp.More(); faceExp.Next() )
5914       {
5915         if ( faceExp.Current().Orientation() != TopAbs_INTERNAL )
5916           continue;
5917
5918         SMESHDS_SubMesh* sm = meshDS->MeshElements( faceExp.Current() );
5919         if ( !sm ) continue;
5920
5921         TIDSortedElemSet facesToOrient;
5922         for ( SMDS_ElemIteratorPtr fIt = sm->GetElements(); fIt->more(); )
5923           facesToOrient.insert( facesToOrient.end(), fIt->next() );
5924         if ( facesToOrient.size() < 2 )
5925           continue;
5926
5927         gp_Dir direction(1,0,0);
5928         TIDSortedElemSet refFaces;
5929         editor.Reorient2D( facesToOrient, direction, refFaces, /*allowNonManifold=*/true );
5930       }
5931     }
5932     return;
5933   }
5934
5935   //================================================================================
5936   /*!
5937    * \brief Create mesh segments.
5938    */
5939   void Hexahedron::addSegments( SMESH_MesherHelper&                      helper,
5940                                 const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap )
5941   {
5942     SMESHDS_Mesh* mesh = helper.GetMeshDS();
5943
5944     std::vector<const SMDS_MeshNode*> nodes;
5945     std::vector<const SMDS_MeshElement *> elems;
5946     map< TGeomID, vector< TGeomID > >::const_iterator e2ff = edge2faceIDsMap.begin();
5947     for ( ; e2ff != edge2faceIDsMap.end(); ++e2ff )
5948     {
5949       const TopoDS_Edge& edge = TopoDS::Edge( _grid->Shape( e2ff->first ));
5950       const TopoDS_Face& face = TopoDS::Face( _grid->Shape( e2ff->second[0] ));
5951       StdMeshers_FaceSide side( face, edge, helper.GetMesh(), /*isFwd=*/true, /*skipMed=*/true );
5952       nodes = side.GetOrderedNodes();
5953
5954       elems.clear();
5955       if ( nodes.size() == 2 )
5956         // check that there is an element connecting two nodes
5957         if ( !mesh->GetElementsByNodes( nodes, elems ))
5958           continue;
5959
5960       for ( size_t i = 1; i < nodes.size(); i++ )
5961       {
5962         if ( mesh->FindEdge( nodes[i-1], nodes[i] ))
5963           continue;
5964         SMDS_MeshElement* segment = mesh->AddEdge( nodes[i-1], nodes[i] );
5965         mesh->SetMeshElementOnShape( segment, e2ff->first );
5966       }
5967     }
5968     return;
5969   }
5970
5971   //================================================================================
5972   /*!
5973    * \brief Return created volumes and volumes that can have free facet because of
5974    *        skipped small volume. Also create mesh faces on free facets
5975    *        of adjacent not-cut volumes if the result volume is too small.
5976    */
5977   void Hexahedron::getBoundaryElems( vector< const SMDS_MeshElement* > & boundaryElems )
5978   {
5979     if ( _hasTooSmall /*|| _volumeDefs.IsEmpty()*/ )
5980     {
5981       // create faces around a missing small volume
5982       TGeomID faceID = 0;
5983       SMESH_MeshEditor editor( _grid->_helper->GetMesh() );
5984       SMESH_MeshEditor::ElemFeatures polygon( SMDSAbs_Face );
5985       SMESHDS_Mesh* meshDS = _grid->_helper->GetMeshDS();
5986       std::vector<const SMDS_MeshElement *> adjVolumes(2);
5987       for ( size_t iF = 0; iF < _polygons.size(); ++iF )
5988       {
5989         const size_t nbLinks = _polygons[ iF ]._links.size();
5990         if ( nbLinks != 4 ) continue;
5991         polygon.myNodes.resize( nbLinks );
5992         polygon.myNodes.back() = 0;
5993         for ( size_t iL = 0, iN = nbLinks - 1; iL < nbLinks; ++iL, --iN )
5994           if ( ! ( polygon.myNodes[iN] = _polygons[ iF ]._links[ iL ].FirstNode()->Node() ))
5995             break;
5996         if ( !polygon.myNodes.back() )
5997           continue;
5998
5999         meshDS->GetElementsByNodes( polygon.myNodes, adjVolumes, SMDSAbs_Volume );
6000         if ( adjVolumes.size() != 1 )
6001           continue;
6002         if ( !adjVolumes[0]->isMarked() )
6003         {
6004           boundaryElems.push_back( adjVolumes[0] );
6005           adjVolumes[0]->setIsMarked( true );
6006         }
6007
6008         bool sameShape = true;
6009         TGeomID shapeID = polygon.myNodes[0]->GetShapeID();
6010         for ( size_t i = 1; i < polygon.myNodes.size() && sameShape; ++i )
6011           sameShape = ( shapeID == polygon.myNodes[i]->GetShapeID() );
6012
6013         if ( !sameShape || !_grid->IsSolid( shapeID ))
6014           continue; // some of shapes must be FACE
6015
6016         if ( !faceID )
6017         {
6018           faceID = getAnyFace();
6019           if ( !faceID )
6020             break;
6021           if ( _grid->IsInternal( faceID ) ||
6022                _grid->IsShared( faceID ) //||
6023                //_grid->IsBoundaryFace( faceID ) -- commented for #19887
6024                ) 
6025             break; // create only if a new face will be used by other 3D algo
6026         }
6027
6028         Solid * solid = _grid->GetOneOfSolids( adjVolumes[0]->GetShapeID() );
6029         if ( !solid->IsOutsideOriented( faceID ))
6030           std::reverse( polygon.myNodes.begin(), polygon.myNodes.end() );
6031
6032         //polygon.SetPoly( polygon.myNodes.size() > 4 );
6033         const SMDS_MeshElement* newFace = editor.AddElement( polygon.myNodes, polygon );
6034         meshDS->SetMeshElementOnShape( newFace, faceID );
6035       }
6036     }
6037
6038     // return created volumes
6039     for ( _volumeDef* volDef = &_volumeDefs; volDef; volDef = volDef->_next )
6040     {
6041       if ( volDef ->_volume &&
6042            !volDef->_volume->IsNull() &&
6043            !volDef->_volume->isMarked() )
6044       {
6045         volDef->_volume->setIsMarked( true );
6046         boundaryElems.push_back( volDef->_volume );
6047
6048         if ( _grid->IsToCheckNodePos() ) // un-mark nodes marked in addVolumes()
6049           for ( size_t iN = 0; iN < volDef->_nodes.size(); ++iN )
6050             volDef->_nodes[iN].Node()->setIsMarked( false );
6051       }
6052       if ( volDef->_brotherVolume.size() > 0 )
6053       {
6054         for (auto _bro : volDef->_brotherVolume )
6055         {
6056           _bro->setIsMarked( true );
6057           boundaryElems.push_back( _bro );
6058         }        
6059       }
6060     }
6061   }
6062
6063   //================================================================================
6064   /*!
6065    * \brief Remove edges and nodes dividing a hexa side in the case if an adjacent
6066    *        volume also sharing the dividing edge is missing due to its small side.
6067    *        Issue #19887.
6068    */
6069   //================================================================================
6070
6071   void Hexahedron::removeExcessSideDivision(const vector< Hexahedron* >& allHexa)
6072   {
6073     if ( ! _volumeDefs.IsPolyhedron() )
6074       return; // not a polyhedron
6075       
6076     // look for a divided side adjacent to a small hexahedron
6077
6078     int di[6] = { 0, 0, 0, 0,-1, 1 };
6079     int dj[6] = { 0, 0,-1, 1, 0, 0 };
6080     int dk[6] = {-1, 1, 0, 0, 0, 0 };
6081
6082     for ( int iF = 0; iF < 6; ++iF ) // loop on 6 sides of a hexahedron
6083     {
6084       size_t neighborIndex = _grid->CellIndex( _i + di[iF],
6085                                                _j + dj[iF],
6086                                                _k + dk[iF] );
6087       if ( neighborIndex >= allHexa.size() ||
6088            !allHexa[ neighborIndex ]       ||
6089            !allHexa[ neighborIndex ]->_hasTooSmall )
6090         continue;
6091
6092       // check if a side is divided into several polygons
6093       for ( _volumeDef* volDef = &_volumeDefs; volDef; volDef = volDef->_next )
6094       {
6095         int nbPolygons = 0, nbNodes = 0;
6096         for ( size_t i = 0; i < volDef->_names.size(); ++i )
6097           if ( volDef->_names[ i ] == _hexQuads[ iF ]._name )
6098           {
6099             ++nbPolygons;
6100             nbNodes += volDef->_quantities[ i ];
6101           }
6102         if ( nbPolygons < 2 )
6103           continue;
6104
6105         // construct loops from polygons
6106         typedef _volumeDef::_linkDef TLinkDef;
6107         std::vector< TLinkDef* > loops;
6108         std::vector< TLinkDef > links( nbNodes );
6109         for ( size_t i = 0, iN = 0, iLoop = 0; iLoop < volDef->_quantities.size(); ++iLoop )
6110         {
6111           size_t nbLinks = volDef->_quantities[ iLoop ];
6112           if ( volDef->_names[ iLoop ] != _hexQuads[ iF ]._name )
6113           {
6114             iN += nbLinks;
6115             continue;
6116           }
6117           loops.push_back( & links[i] );
6118           for ( size_t n = 0; n < nbLinks-1; ++n, ++i, ++iN )
6119           {
6120             links[i].init( volDef->_nodes[iN], volDef->_nodes[iN+1], iLoop );
6121             links[i].setNext( &links[i+1] );
6122           }
6123           links[i].init( volDef->_nodes[iN], volDef->_nodes[iN-nbLinks+1], iLoop );
6124           links[i].setNext( &links[i-nbLinks+1] );
6125           ++i; ++iN;
6126         }
6127
6128         // look for equal links in different loops and join such loops
6129         bool loopsJoined = false;
6130         std::set< TLinkDef > linkSet;
6131         for ( size_t iLoop = 0; iLoop < loops.size(); ++iLoop )
6132         {
6133           TLinkDef* beg = 0;
6134           for ( TLinkDef* l = loops[ iLoop ]; l != beg; l = l->_next ) // walk around the iLoop
6135           {
6136             std::pair< std::set< TLinkDef >::iterator, bool > it2new = linkSet.insert( *l );
6137             if ( !it2new.second ) // equal found, join loops
6138             {
6139               const TLinkDef* equal = &(*it2new.first);
6140               if ( equal->_loopIndex == l->_loopIndex )
6141                 continue; // error?
6142
6143               loopsJoined = true;
6144
6145               for ( size_t i = iLoop - 1; i < loops.size(); --i )
6146                 if ( loops[ i ] && loops[ i ]->_loopIndex == equal->_loopIndex )
6147                   loops[ i ] = 0;
6148
6149               // exclude l and equal and join two loops
6150               if ( l->_prev != equal )
6151                 l->_prev->setNext( equal->_next );
6152               if ( equal->_prev != l )
6153                 equal->_prev->setNext( l->_next );
6154
6155               if ( volDef->_quantities[ l->_loopIndex ] > 0 )
6156                 volDef->_quantities[ l->_loopIndex     ] *= -1;
6157               if ( volDef->_quantities[ equal->_loopIndex ] > 0 )
6158                 volDef->_quantities[ equal->_loopIndex ] *= -1;
6159
6160               if ( loops[ iLoop ] == l )
6161                 loops[ iLoop ] = l->_prev->_next;
6162             }
6163             beg = loops[ iLoop ];
6164           }
6165         }
6166         // update volDef
6167         if ( loopsJoined )
6168         {
6169           // set unchanged polygons
6170           std::vector< int >             newQuantities;
6171           std::vector< _volumeDef::_nodeDef > newNodes;
6172           vector< SMESH_Block::TShapeID >     newNames;
6173           newQuantities.reserve( volDef->_quantities.size() );
6174           newNodes.reserve     ( volDef->_nodes.size() );
6175           newNames.reserve     ( volDef->_names.size() );
6176           for ( size_t i = 0, iLoop = 0; iLoop < volDef->_quantities.size(); ++iLoop )
6177           {
6178             if ( volDef->_quantities[ iLoop ] < 0 )
6179             {
6180               i -= volDef->_quantities[ iLoop ];
6181               continue;
6182             }
6183             newQuantities.push_back( volDef->_quantities[ iLoop ]);
6184             newNodes.insert( newNodes.end(),
6185                              volDef->_nodes.begin() + i,
6186                              volDef->_nodes.begin() + i + newQuantities.back() );
6187             newNames.push_back( volDef->_names[ iLoop ]);
6188             i += volDef->_quantities[ iLoop ];
6189           }
6190
6191           // set joined loops
6192           for ( size_t iLoop = 0; iLoop < loops.size(); ++iLoop )
6193           {
6194             if ( !loops[ iLoop ] )
6195               continue;
6196             newQuantities.push_back( 0 );
6197             TLinkDef* beg = 0;
6198             for ( TLinkDef* l = loops[ iLoop ]; l != beg; l = l->_next, ++newQuantities.back() )
6199             {
6200               newNodes.push_back( l->_node1 );
6201               beg = loops[ iLoop ];
6202             }
6203             newNames.push_back( _hexQuads[ iF ]._name );
6204           }
6205           volDef->_quantities.swap( newQuantities );
6206           volDef->_nodes.swap( newNodes );
6207           volDef->_names.swap( newNames );
6208         }
6209       } // loop on volDef's
6210     } // loop on hex sides
6211
6212     return;
6213   } // removeExcessSideDivision()
6214
6215
6216   //================================================================================
6217   /*!
6218    * \brief Remove nodes splitting Cartesian cell edges in the case if a node
6219    *        is used in every cells only by two polygons sharing the edge
6220    *        Issue #19887.
6221    */
6222   //================================================================================
6223
6224   void Hexahedron::removeExcessNodes(vector< Hexahedron* >& allHexa)
6225   {
6226     if ( ! _volumeDefs.IsPolyhedron() )
6227       return; // not a polyhedron
6228
6229     typedef vector< _volumeDef::_nodeDef >::iterator TNodeIt;
6230     vector< int > nodesInPoly[ 4 ]; // node index in _volumeDefs._nodes
6231     vector< int > volDefInd  [ 4 ]; // index of a _volumeDefs
6232     Hexahedron*   hexa       [ 4 ];
6233     int i,j,k, cellIndex, iLink = 0, iCellLink;
6234     for ( int iDir = 0; iDir < 3; ++iDir )
6235     {
6236       CellsAroundLink fourCells( _grid, iDir );
6237       for ( int iL = 0; iL < 4; ++iL, ++iLink ) // 4 links in a direction
6238       {
6239         _Link& link = _hexLinks[ iLink ];
6240         fourCells.Init( _i, _j, _k, iLink );
6241
6242         for ( size_t iP = 0; iP < link._fIntPoints.size(); ++iP ) // loop on nodes on the link
6243         {
6244           bool nodeRemoved = true;
6245           _volumeDef::_nodeDef node; node._intPoint = link._fIntPoints[iP];
6246
6247           for ( size_t i = 0, nb = _volumeDefs.size(); i < nb &&  nodeRemoved; ++i )
6248             if ( _volumeDef* vol = _volumeDefs.at( i ))
6249               nodeRemoved =
6250                 ( std::find( vol->_nodes.begin(), vol->_nodes.end(), node ) == vol->_nodes.end() );
6251           if ( nodeRemoved )
6252             continue; // node already removed
6253
6254           // check if a node encounters zero or two times in 4 cells sharing iLink
6255           // if so, the node can be removed from the cells
6256           bool       nodeIsOnEdge = true;
6257           int nbPolyhedraWithNode = 0;
6258           for ( int iC = 0; iC < 4; ++iC ) // loop on 4 cells sharing a link
6259           {
6260             nodesInPoly[ iC ].clear();
6261             volDefInd  [ iC ].clear();
6262             hexa       [ iC ] = 0;
6263             if ( !fourCells.GetCell( iC, i,j,k, cellIndex, iCellLink ))
6264               continue;
6265             hexa[ iC ] = allHexa[ cellIndex ];
6266             if ( !hexa[ iC ])
6267               continue;
6268             for ( size_t i = 0, nb = hexa[ iC ]->_volumeDefs.size(); i < nb; ++i )
6269               if ( _volumeDef* vol = hexa[ iC ]->_volumeDefs.at( i ))
6270               {
6271                 for ( TNodeIt nIt = vol->_nodes.begin(); nIt != vol->_nodes.end(); ++nIt )
6272                 {
6273                   nIt = std::find( nIt, vol->_nodes.end(), node );
6274                   if ( nIt != vol->_nodes.end() )
6275                   {
6276                     nodesInPoly[ iC ].push_back( std::distance( vol->_nodes.begin(), nIt ));
6277                     volDefInd  [ iC ].push_back( i );
6278                   }
6279                   else
6280                     break;
6281                 }
6282                 nbPolyhedraWithNode += ( !nodesInPoly[ iC ].empty() );
6283               }
6284             if ( nodesInPoly[ iC ].size() != 0 &&
6285                  nodesInPoly[ iC ].size() != 2 )
6286             {
6287               nodeIsOnEdge = false;
6288               break;
6289             }
6290           } // loop  on 4 cells
6291
6292           // remove nodes from polyhedra
6293           if ( nbPolyhedraWithNode > 0 && nodeIsOnEdge )
6294           {
6295             for ( int iC = 0; iC < 4; ++iC ) // loop on 4 cells sharing the link
6296             {
6297               if ( nodesInPoly[ iC ].empty() )
6298                 continue;
6299               for ( int i = volDefInd[ iC ].size() - 1; i >= 0; --i )
6300               {
6301                 _volumeDef* vol = hexa[ iC ]->_volumeDefs.at( volDefInd[ iC ][ i ]);
6302                 int nIndex = nodesInPoly[ iC ][ i ];
6303                 // decrement _quantities
6304                 for ( size_t iQ = 0; iQ < vol->_quantities.size(); ++iQ )
6305                   if ( nIndex < vol->_quantities[ iQ ])
6306                   {
6307                     vol->_quantities[ iQ ]--;
6308                     break;
6309                   }
6310                   else
6311                   {
6312                     nIndex -= vol->_quantities[ iQ ];
6313                   }
6314                 vol->_nodes.erase( vol->_nodes.begin() + nodesInPoly[ iC ][ i ]);
6315
6316                 if ( i == 0 &&
6317                      vol->_nodes.size() == 6 * 4 &&
6318                      vol->_quantities.size() == 6 ) // polyhedron becomes hexahedron?
6319                 {
6320                   bool allQuads = true;
6321                   for ( size_t iQ = 0; iQ < vol->_quantities.size() &&  allQuads; ++iQ )
6322                     allQuads = ( vol->_quantities[ iQ ] == 4 );
6323                   if ( allQuads )
6324                   {
6325                     // set side nodes as this: bottom, top, top, ...
6326                     int iTop = 0, iBot = 0; // side indices
6327                     for ( int iS = 0; iS < 6; ++iS )
6328                     {
6329                       if ( vol->_names[ iS ] == SMESH_Block::ID_Fxy0 )
6330                         iBot = iS;
6331                       if ( vol->_names[ iS ] == SMESH_Block::ID_Fxy1 )
6332                         iTop = iS;
6333                     }
6334                     if ( iBot != 0 )
6335                     {
6336                       if ( iTop == 0 )
6337                       {
6338                         std::copy( vol->_nodes.begin(),
6339                                    vol->_nodes.begin() + 4,
6340                                    vol->_nodes.begin() + 4 );
6341                         iTop = 1;
6342                       }
6343                       std::copy( vol->_nodes.begin() + 4 * iBot,
6344                                  vol->_nodes.begin() + 4 * ( iBot + 1),
6345                                  vol->_nodes.begin() );
6346                     }
6347                     if ( iTop != 1 )
6348                       std::copy( vol->_nodes.begin() + 4 * iTop,
6349                                  vol->_nodes.begin() + 4 * ( iTop + 1),
6350                                  vol->_nodes.begin() + 4 );
6351
6352                     std::copy( vol->_nodes.begin() + 4,
6353                                vol->_nodes.begin() + 8,
6354                                vol->_nodes.begin() + 8 );
6355                     // set up top facet nodes by comparing their uvw with bottom nodes
6356                     E_IntersectPoint ip[8];
6357                     for ( int iN = 0; iN < 8; ++iN )
6358                     {
6359                       SMESH_NodeXYZ p = vol->_nodes[ iN ].Node();
6360                       _grid->ComputeUVW( p, ip[ iN ]._uvw );
6361                     }
6362                     const double tol2 = _grid->_tol * _grid->_tol;
6363                     for ( int iN = 0; iN < 4; ++iN )
6364                     {
6365                       gp_Pnt2d pBot( ip[ iN ]._uvw[0], ip[ iN ]._uvw[1] );
6366                       for ( int iT = 4; iT < 8; ++iT )
6367                       {
6368                         gp_Pnt2d pTop( ip[ iT ]._uvw[0], ip[ iT ]._uvw[1] );
6369                         if ( pBot.SquareDistance( pTop ) < tol2 )
6370                         {
6371                           // vol->_nodes[ iN + 4 ]._node = ip[ iT ]._node;
6372                           // vol->_nodes[ iN + 4 ]._intPoint = 0;
6373                           vol->_nodes[ iN + 4 ] = vol->_nodes[ iT + 4 ];
6374                           break;
6375                         }
6376                       }
6377                     }
6378                     vol->_nodes.resize( 8 );
6379                     vol->_quantities.clear();
6380                     //vol->_names.clear();
6381                   }
6382                 }
6383               } // loop on _volumeDefs
6384             } // loop on 4 cell abound a link
6385           } // if ( nodeIsOnEdge )
6386         } // loop on intersection points of a link
6387       } // loop on 4 links of a direction
6388     } // loop on 3 directions
6389
6390     return;
6391
6392   } // removeExcessNodes()
6393
6394   //================================================================================
6395   /*!
6396    * \brief [Issue #19913] Modify _hexLinks._splits to prevent creating overlapping volumes
6397    */
6398   //================================================================================
6399
6400   void Hexahedron::preventVolumesOverlapping()
6401   {
6402     // Cut off a quadrangle corner if two links sharing the corner
6403     // are shared by same two solids, in this case each of solids gets
6404     // a triangle for it-self.
6405     std::vector< TGeomID > soIDs[4];
6406     for ( int iF = 0; iF < 6; ++iF ) // loop on 6 sides of a hexahedron
6407     {
6408       _Face& quad = _hexQuads[ iF ] ;
6409
6410       int iFOpposite = iF + ( iF % 2 ? -1 : 1 );
6411       _Face& quadOpp = _hexQuads[ iFOpposite ] ;
6412
6413       int nbSides = 0, nbSidesOpp = 0;
6414       for ( int iE = 0; iE < 4; ++iE ) // loop on 4 sides of a quadrangle
6415       {
6416         nbSides    += ( quad._links   [ iE ].NbResultLinks() > 0 );
6417         nbSidesOpp += ( quadOpp._links[ iE ].NbResultLinks() > 0 );
6418       }
6419       if ( nbSides < 4 || nbSidesOpp != 2 )
6420         continue;
6421
6422       for ( int iE = 0; iE < 4; ++iE )
6423       {
6424         soIDs[ iE ].clear();
6425         _Node* n = quad._links[ iE ].FirstNode();
6426         if ( n->_intPoint && n->_intPoint->_faceIDs.size() )
6427           soIDs[ iE ] = _grid->GetSolidIDs( n->_intPoint->_faceIDs[0] );
6428       }
6429       if ((( soIDs[0].size() >= 2 ) +
6430            ( soIDs[1].size() >= 2 ) +
6431            ( soIDs[2].size() >= 2 ) +
6432            ( soIDs[3].size() >= 2 ) ) < 3 )
6433         continue;
6434
6435       bool done = false;
6436       for ( int i = 0; i < 4; ++i )
6437       {
6438         int i1 = _grid->_helper->WrapIndex( i + 1, 4 );
6439         int i2 = _grid->_helper->WrapIndex( i + 2, 4 );
6440         int i3 = _grid->_helper->WrapIndex( i + 3, 4 );
6441         if ( soIDs[i1].size() == 2 && soIDs[i ] != soIDs[i1] &&
6442              soIDs[i2].size() == 2 && soIDs[i1] == soIDs[i2] &&
6443              soIDs[i3].size() == 2 && soIDs[i2] == soIDs[i3] )
6444         {
6445           quad._links[ i1 ]._link->_splits.clear();
6446           quad._links[ i2 ]._link->_splits.clear();
6447           done = true;
6448           break;
6449         }
6450       }
6451       if ( done )
6452         break;
6453     }
6454     return;
6455   } // preventVolumesOverlapping()
6456
6457   //================================================================================
6458   /*!
6459    * \brief Set to _hexLinks a next portion of splits located on one side of INTERNAL FACEs
6460    */
6461   bool Hexahedron::_SplitIterator::Next()
6462   {
6463     if ( _iterationNb > 0 )
6464       // count used splits
6465       for ( size_t i = 0; i < _splits.size(); ++i )
6466       {
6467         if ( _splits[i]._iCheckIteration == _iterationNb )
6468         {
6469           _splits[i]._isUsed = _splits[i]._checkedSplit->_faces[1];
6470           _nbUsed += _splits[i]._isUsed;
6471         }
6472         if ( !More() )
6473           return false;
6474       }
6475
6476     ++_iterationNb;
6477
6478     bool toTestUsed = ( _nbChecked >= _splits.size() );
6479     if ( toTestUsed )
6480     {
6481       // all splits are checked; find all not used splits
6482       for ( size_t i = 0; i < _splits.size(); ++i )
6483         if ( !_splits[i].IsCheckedOrUsed( toTestUsed ))
6484           _splits[i]._iCheckIteration = _iterationNb;
6485
6486       _nbUsed = _splits.size(); // to stop iteration
6487     }
6488     else
6489     {
6490       // get any not used/checked split to start from
6491       _freeNodes.clear();
6492       for ( size_t i = 0; i < _splits.size(); ++i )
6493       {
6494         if ( !_splits[i].IsCheckedOrUsed( toTestUsed ))
6495         {
6496           _freeNodes.push_back( _splits[i]._nodes[0] );
6497           _freeNodes.push_back( _splits[i]._nodes[1] );
6498           _splits[i]._iCheckIteration = _iterationNb;
6499           break;
6500         }
6501       }
6502       // find splits connected to the start one via _freeNodes
6503       for ( size_t iN = 0; iN < _freeNodes.size(); ++iN )
6504       {
6505         for ( size_t iS = 0; iS < _splits.size(); ++iS )
6506         {
6507           if ( _splits[iS].IsCheckedOrUsed( toTestUsed ))
6508             continue;
6509           int iN2 = -1;
6510           if (      _freeNodes[iN] == _splits[iS]._nodes[0] )
6511             iN2 = 1;
6512           else if ( _freeNodes[iN] == _splits[iS]._nodes[1] )
6513             iN2 = 0;
6514           else
6515             continue;
6516           if ( _freeNodes[iN]->_isInternalFlags > 0 )
6517           {
6518             if ( _splits[iS]._nodes[ iN2 ]->_isInternalFlags == 0 )
6519               continue;
6520             if ( !_splits[iS]._nodes[ iN2 ]->IsLinked( _freeNodes[iN]->_intPoint ))
6521               continue;
6522           }
6523           _splits[iS]._iCheckIteration = _iterationNb;
6524           _freeNodes.push_back( _splits[iS]._nodes[ iN2 ]);
6525         }
6526       }
6527     }
6528     // set splits to hex links
6529
6530     for ( int iL = 0; iL < 12; ++iL )
6531       _hexLinks[ iL ]._splits.clear();
6532
6533     _Link split;
6534     for ( size_t i = 0; i < _splits.size(); ++i )
6535     {
6536       if ( _splits[i]._iCheckIteration == _iterationNb )
6537       {
6538         split._nodes[0] = _splits[i]._nodes[0];
6539         split._nodes[1] = _splits[i]._nodes[1];
6540         _Link & hexLink = _hexLinks[ _splits[i]._linkID ];
6541         hexLink._splits.push_back( split );
6542         _splits[i]._checkedSplit = & hexLink._splits.back();
6543         ++_nbChecked;
6544       }
6545     }
6546     return More();
6547   }
6548
6549   //================================================================================
6550   /*!
6551    * \brief computes exact bounding box with axes parallel to given ones
6552    */
6553   //================================================================================
6554
6555   void getExactBndBox( const vector< TopoDS_Shape >& faceVec,
6556                        const double*                 axesDirs,
6557                        Bnd_Box&                      shapeBox )
6558   {
6559     BRep_Builder b;
6560     TopoDS_Compound allFacesComp;
6561     b.MakeCompound( allFacesComp );
6562     for ( size_t iF = 0; iF < faceVec.size(); ++iF )
6563       b.Add( allFacesComp, faceVec[ iF ] );
6564
6565     double sP[6]; // aXmin, aYmin, aZmin, aXmax, aYmax, aZmax
6566     shapeBox.Get(sP[0],sP[1],sP[2],sP[3],sP[4],sP[5]);
6567     double farDist = 0;
6568     for ( int i = 0; i < 6; ++i )
6569       farDist = Max( farDist, 10 * sP[i] );
6570
6571     gp_XYZ axis[3] = { gp_XYZ( axesDirs[0], axesDirs[1], axesDirs[2] ),
6572                        gp_XYZ( axesDirs[3], axesDirs[4], axesDirs[5] ),
6573                        gp_XYZ( axesDirs[6], axesDirs[7], axesDirs[8] ) };
6574     axis[0].Normalize();
6575     axis[1].Normalize();
6576     axis[2].Normalize();
6577
6578     gp_Mat basis( axis[0], axis[1], axis[2] );
6579     gp_Mat bi = basis.Inverted();
6580
6581     gp_Pnt pMin, pMax;
6582     for ( int iDir = 0; iDir < 3; ++iDir )
6583     {
6584       gp_XYZ axis0 = axis[ iDir ];
6585       gp_XYZ axis1 = axis[ ( iDir + 1 ) % 3 ];
6586       gp_XYZ axis2 = axis[ ( iDir + 2 ) % 3 ];
6587       for ( int isMax = 0; isMax < 2; ++isMax )
6588       {
6589         double shift = isMax ? farDist : -farDist;
6590         gp_XYZ orig = shift * axis0;
6591         gp_XYZ norm = axis1 ^ axis2;
6592         gp_Pln pln( orig, norm );
6593         norm = pln.Axis().Direction().XYZ();
6594         BRepBuilderAPI_MakeFace plane( pln, -farDist, farDist, -farDist, farDist );
6595
6596         gp_Pnt& pAxis = isMax ? pMax : pMin;
6597         gp_Pnt pPlane, pFaces;
6598         double dist = GEOMUtils::GetMinDistance( plane, allFacesComp, pPlane, pFaces );
6599         if ( dist < 0 )
6600         {
6601           Bnd_B3d bb;
6602           gp_XYZ corner;
6603           for ( int i = 0; i < 2; ++i ) {
6604             corner.SetCoord( 1, sP[ i*3 ]);
6605             for ( int j = 0; j < 2; ++j ) {
6606               corner.SetCoord( 2, sP[ i*3 + 1 ]);
6607               for ( int k = 0; k < 2; ++k )
6608               {
6609                 corner.SetCoord( 3, sP[ i*3 + 2 ]);
6610                 corner *= bi;
6611                 bb.Add( corner );
6612               }
6613             }
6614           }
6615           corner = isMax ? bb.CornerMax() : bb.CornerMin();
6616           pAxis.SetCoord( iDir+1, corner.Coord( iDir+1 ));
6617         }
6618         else
6619         {
6620           gp_XYZ pf = pFaces.XYZ() * bi;
6621           pAxis.SetCoord( iDir+1, pf.Coord( iDir+1 ) );
6622         }
6623       }
6624     } // loop on 3 axes
6625
6626     shapeBox.SetVoid();
6627     shapeBox.Add( pMin );
6628     shapeBox.Add( pMax );
6629
6630     return;
6631   }
6632
6633 } // namespace
6634
6635 //=============================================================================
6636 /*!
6637  * \brief Generates 3D structured Cartesian mesh in the internal part of
6638  * solid shapes and polyhedral volumes near the shape boundary.
6639  *  \param theMesh - mesh to fill in
6640  *  \param theShape - a compound of all SOLIDs to mesh
6641  *  \retval bool - true in case of success
6642  */
6643 //=============================================================================
6644
6645 bool StdMeshers_Cartesian_3D::Compute(SMESH_Mesh &         theMesh,
6646                                       const TopoDS_Shape & theShape)
6647 {
6648   if ( _hypViscousLayers )
6649   {
6650     const StdMeshers_ViscousLayers* hypViscousLayers = _hypViscousLayers;
6651     _hypViscousLayers = nullptr;
6652
6653     StdMeshers_Cartesian_VL::ViscousBuilder builder( hypViscousLayers, theMesh, theShape );
6654
6655     std::string error;
6656     TopoDS_Shape offsetShape = builder.MakeOffsetShape( theShape, theMesh, error );
6657     if ( offsetShape.IsNull() )
6658       throw SALOME_Exception( error );
6659
6660     SMESH_Mesh* offsetMesh = new TmpMesh(); 
6661     offsetMesh->ShapeToMesh( offsetShape );
6662     offsetMesh->GetSubMesh( offsetShape )->DependsOn();
6663
6664     this->_isComputeOffset = true;
6665     if ( ! this->Compute( *offsetMesh, offsetShape ))
6666       return false;
6667
6668     return builder.MakeViscousLayers( *offsetMesh, theMesh, theShape );
6669   }
6670
6671   // The algorithm generates the mesh in following steps:
6672
6673   // 1) Intersection of grid lines with the geometry boundary.
6674   // This step allows to find out if a given node of the initial grid is
6675   // inside or outside the geometry.
6676
6677   // 2) For each cell of the grid, check how many of it's nodes are outside
6678   // of the geometry boundary. Depending on a result of this check
6679   // - skip a cell, if all it's nodes are outside
6680   // - skip a cell, if it is too small according to the size threshold
6681   // - add a hexahedron in the mesh, if all nodes are inside
6682   // - add a polyhedron in the mesh, if some nodes are inside and some outside
6683
6684   _computeCanceled = false;
6685
6686   SMESH_MesherHelper helper( theMesh );
6687   SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
6688
6689   try
6690   {
6691     Grid grid;
6692     grid._helper                         = &helper;
6693     grid._toAddEdges                     = _hyp->GetToAddEdges();
6694     grid._toCreateFaces                  = _hyp->GetToCreateFaces();
6695     grid._toConsiderInternalFaces        = _hyp->GetToConsiderInternalFaces();
6696     grid._toUseThresholdForInternalFaces = _hyp->GetToUseThresholdForInternalFaces();
6697     grid._sizeThreshold                  = _hyp->GetSizeThreshold();
6698     grid._toUseQuanta                    = _hyp->GetToUseQuanta();
6699     grid._quanta                         = _hyp->GetQuanta();
6700     if ( _isComputeOffset )
6701     {
6702       grid._toAddEdges = true;
6703       grid._toCreateFaces = true;
6704     }
6705     grid.InitGeometry( theShape );
6706
6707     vector< TopoDS_Shape > faceVec;
6708     {
6709       TopTools_MapOfShape faceMap;
6710       TopExp_Explorer fExp;
6711       for ( fExp.Init( theShape, TopAbs_FACE ); fExp.More(); fExp.Next() )
6712       {
6713         bool isNewFace = faceMap.Add( fExp.Current() );
6714         if ( !grid._toConsiderInternalFaces )
6715           if ( !isNewFace || fExp.Current().Orientation() == TopAbs_INTERNAL )
6716             // remove an internal face
6717             faceMap.Remove( fExp.Current() );
6718       }
6719       faceVec.reserve( faceMap.Extent() );
6720       faceVec.assign( faceMap.cbegin(), faceMap.cend() );
6721     }
6722     vector<FaceGridIntersector> facesItersectors( faceVec.size() );
6723     Bnd_Box shapeBox;
6724     for ( size_t i = 0; i < faceVec.size(); ++i )
6725     {
6726       facesItersectors[i]._face   = TopoDS::Face( faceVec[i] );
6727       facesItersectors[i]._faceID = grid.ShapeID( faceVec[i] );
6728       facesItersectors[i]._grid   = &grid;
6729       shapeBox.Add( facesItersectors[i].GetFaceBndBox() );
6730     }
6731     getExactBndBox( faceVec, _hyp->GetAxisDirs(), shapeBox );
6732
6733
6734     vector<double> xCoords, yCoords, zCoords;
6735     _hyp->GetCoordinates( xCoords, yCoords, zCoords, shapeBox );
6736
6737     grid.SetCoordinates( xCoords, yCoords, zCoords, _hyp->GetAxisDirs(), shapeBox );
6738
6739     if ( _computeCanceled ) return false;
6740
6741 #ifdef WITH_TBB
6742     { // copy partner faces and curves of not thread-safe types
6743       set< const Standard_Transient* > tshapes;
6744       BRepBuilderAPI_Copy copier;
6745       for ( size_t i = 0; i < facesItersectors.size(); ++i )
6746       {
6747         if ( !facesItersectors[i].IsThreadSafe( tshapes ))
6748         {
6749           copier.Perform( facesItersectors[i]._face );
6750           facesItersectors[i]._face = TopoDS::Face( copier );
6751         }
6752       }
6753     }
6754     // Intersection of grid lines with the geometry boundary.
6755     tbb::parallel_for ( tbb::blocked_range<size_t>( 0, facesItersectors.size() ),
6756                         ParallelIntersector( facesItersectors ),
6757                         tbb::simple_partitioner());
6758 #else
6759     for ( size_t i = 0; i < facesItersectors.size(); ++i )
6760       facesItersectors[i].Intersect();
6761 #endif
6762
6763     // put intersection points onto the GridLine's; this is done after intersection
6764     // to avoid contention of facesItersectors for writing into the same GridLine
6765     // in case of parallel work of facesItersectors
6766     for ( size_t i = 0; i < facesItersectors.size(); ++i )
6767       facesItersectors[i].StoreIntersections();
6768
6769     if ( _computeCanceled ) return false;
6770
6771     // create nodes on the geometry
6772     grid.ComputeNodes( helper );
6773
6774     if ( _computeCanceled ) return false;
6775
6776     // get EDGEs to take into account
6777     map< TGeomID, vector< TGeomID > > edge2faceIDsMap;
6778     grid.GetEdgesToImplement( edge2faceIDsMap, theShape, faceVec );
6779
6780     // create volume elements
6781     Hexahedron hex( &grid );
6782     int nbAdded = hex.MakeElements( helper, edge2faceIDsMap );
6783
6784     if ( nbAdded > 0 )
6785     {
6786       if ( !grid._toConsiderInternalFaces )
6787       {
6788         // make all SOLIDs computed
6789         TopExp_Explorer solidExp( theShape, TopAbs_SOLID );
6790         if ( SMESHDS_SubMesh* sm1 = meshDS->MeshElements( solidExp.Current()) )
6791         {
6792           SMDS_ElemIteratorPtr volIt = sm1->GetElements();
6793           for ( ; solidExp.More() && volIt->more(); solidExp.Next() )
6794           {
6795             const SMDS_MeshElement* vol = volIt->next();
6796             sm1->RemoveElement( vol );
6797             meshDS->SetMeshElementOnShape( vol, solidExp.Current() );
6798           }
6799         }
6800       }
6801       // make other sub-shapes computed
6802       setSubmeshesComputed( theMesh, theShape );
6803     }
6804
6805     // remove free nodes
6806     //if ( SMESHDS_SubMesh * smDS = meshDS->MeshElements( helper.GetSubShapeID() ))
6807     {
6808       std::vector< const SMDS_MeshNode* > nodesToRemove;
6809       // get intersection nodes
6810       for ( int iDir = 0; iDir < 3; ++iDir )
6811       {
6812         vector< GridLine >& lines = grid._lines[ iDir ];
6813         for ( size_t i = 0; i < lines.size(); ++i )
6814         {
6815           multiset< F_IntersectPoint >::iterator ip = lines[i]._intPoints.begin();
6816           for ( ; ip != lines[i]._intPoints.end(); ++ip )
6817             if ( ip->_node &&
6818                  !ip->_node->IsNull() &&
6819                  ip->_node->NbInverseElements() == 0 &&
6820                  !ip->_node->isMarked() )
6821             {
6822               nodesToRemove.push_back( ip->_node );
6823               ip->_node->setIsMarked( true );
6824             }
6825         }
6826       }
6827       // get grid nodes
6828       for ( size_t i = 0; i < grid._nodes.size(); ++i )
6829         if ( grid._nodes[i] &&
6830              !grid._nodes[i]->IsNull() &&
6831              grid._nodes[i]->NbInverseElements() == 0 &&
6832              !grid._nodes[i]->isMarked() )
6833         {
6834           nodesToRemove.push_back( grid._nodes[i] );
6835           grid._nodes[i]->setIsMarked( true );
6836         }
6837
6838       for ( size_t i = 0; i < grid._allBorderNodes.size(); ++i )
6839         if ( grid._allBorderNodes[i] &&
6840              !grid._allBorderNodes[i]->IsNull() &&
6841              grid._allBorderNodes[i]->NbInverseElements() == 0 )
6842         {
6843           nodesToRemove.push_back( grid._allBorderNodes[i] );
6844           grid._allBorderNodes[i]->setIsMarked( true );
6845         }
6846
6847       // do remove
6848       for ( size_t i = 0; i < nodesToRemove.size(); ++i )
6849         meshDS->RemoveFreeNode( nodesToRemove[i], /*smD=*/0, /*fromGroups=*/false );
6850     }
6851
6852     return nbAdded;
6853
6854   }
6855   // SMESH_ComputeError is not caught at SMESH_submesh level for an unknown reason
6856   catch ( SMESH_ComputeError& e)
6857   {
6858     return error( SMESH_ComputeErrorPtr( new SMESH_ComputeError( e )));
6859   }
6860   return false;
6861 }
6862
6863 //=============================================================================
6864 /*!
6865  *  Evaluate
6866  */
6867 //=============================================================================
6868
6869 bool StdMeshers_Cartesian_3D::Evaluate(SMESH_Mesh &         /*theMesh*/,
6870                                        const TopoDS_Shape & /*theShape*/,
6871                                        MapShapeNbElems&     /*theResMap*/)
6872 {
6873   // TODO
6874 //   std::vector<int> aResVec(SMDSEntity_Last);
6875 //   for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
6876 //   if(IsQuadratic) {
6877 //     aResVec[SMDSEntity_Quad_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
6878 //     int nb1d_face0_int = ( nb2d_face0*4 - nb1d ) / 2;
6879 //     aResVec[SMDSEntity_Node] = nb0d_face0 * ( 2*nb2d/nb1d - 1 ) - nb1d_face0_int * nb2d/nb1d;
6880 //   }
6881 //   else {
6882 //     aResVec[SMDSEntity_Node] = nb0d_face0 * ( nb2d/nb1d - 1 );
6883 //     aResVec[SMDSEntity_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
6884 //   }
6885 //   SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
6886 //   aResMap.insert(std::make_pair(sm,aResVec));
6887
6888   return true;
6889 }
6890
6891 //=============================================================================
6892 namespace
6893 {
6894   /*!
6895    * \brief Event listener setting/unsetting _alwaysComputed flag to
6896    *        submeshes of inferior levels to prevent their computing
6897    */
6898   struct _EventListener : public SMESH_subMeshEventListener
6899   {
6900     string _algoName;
6901
6902     _EventListener(const string& algoName):
6903       SMESH_subMeshEventListener(/*isDeletable=*/true,"StdMeshers_Cartesian_3D::_EventListener"),
6904       _algoName(algoName)
6905     {}
6906     // --------------------------------------------------------------------------------
6907     // setting/unsetting _alwaysComputed flag to submeshes of inferior levels
6908     //
6909     static void setAlwaysComputed( const bool     isComputed,
6910                                    SMESH_subMesh* subMeshOfSolid)
6911     {
6912       SMESH_subMeshIteratorPtr smIt =
6913         subMeshOfSolid->getDependsOnIterator(/*includeSelf=*/false, /*complexShapeFirst=*/false);
6914       while ( smIt->more() )
6915       {
6916         SMESH_subMesh* sm = smIt->next();
6917         sm->SetIsAlwaysComputed( isComputed );
6918       }
6919       subMeshOfSolid->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
6920     }
6921
6922     // --------------------------------------------------------------------------------
6923     // unsetting _alwaysComputed flag if "Cartesian_3D" was removed
6924     //
6925     virtual void ProcessEvent(const int          /*event*/,
6926                               const int          eventType,
6927                               SMESH_subMesh*     subMeshOfSolid,
6928                               SMESH_subMeshEventListenerData* /*data*/,
6929                               const SMESH_Hypothesis*         /*hyp*/ = 0)
6930     {
6931       if ( eventType == SMESH_subMesh::COMPUTE_EVENT )
6932       {
6933         setAlwaysComputed( subMeshOfSolid->GetComputeState() == SMESH_subMesh::COMPUTE_OK,
6934                            subMeshOfSolid );
6935       }
6936       else
6937       {
6938         SMESH_Algo* algo3D = subMeshOfSolid->GetAlgo();
6939         if ( !algo3D || _algoName != algo3D->GetName() )
6940           setAlwaysComputed( false, subMeshOfSolid );
6941       }
6942     }
6943
6944     // --------------------------------------------------------------------------------
6945     // set the event listener
6946     //
6947     static void SetOn( SMESH_subMesh* subMeshOfSolid, const string& algoName )
6948     {
6949       subMeshOfSolid->SetEventListener( new _EventListener( algoName ),
6950                                         /*data=*/0,
6951                                         subMeshOfSolid );
6952     }
6953
6954   }; // struct _EventListener
6955
6956 } // namespace
6957
6958 //================================================================================
6959 /*!
6960  * \brief Sets event listener to submeshes if necessary
6961  *  \param subMesh - submesh where algo is set
6962  * This method is called when a submesh gets HYP_OK algo_state.
6963  * After being set, event listener is notified on each event of a submesh.
6964  */
6965 //================================================================================
6966
6967 void StdMeshers_Cartesian_3D::SetEventListener(SMESH_subMesh* subMesh)
6968 {
6969   _EventListener::SetOn( subMesh, GetName() );
6970 }
6971
6972 //================================================================================
6973 /*!
6974  * \brief Set _alwaysComputed flag to submeshes of inferior levels to avoid their computing
6975  */
6976 //================================================================================
6977
6978 void StdMeshers_Cartesian_3D::setSubmeshesComputed(SMESH_Mesh&         theMesh,
6979                                                    const TopoDS_Shape& theShape)
6980 {
6981   for ( TopExp_Explorer soExp( theShape, TopAbs_SOLID ); soExp.More(); soExp.Next() )
6982     _EventListener::setAlwaysComputed( true, theMesh.GetSubMesh( soExp.Current() ));
6983 }