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