Salome HOME
22536: EDF 2876 SMESH : Problem with BodyFitting
[modules/smesh.git] / src / StdMeshers / StdMeshers_Cartesian_3D.cxx
1 // Copyright (C) 2007-2014  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
27 #include "SMDS_MeshNode.hxx"
28 #include "SMESH_Block.hxx"
29 #include "SMESH_Comment.hxx"
30 #include "SMESH_Mesh.hxx"
31 #include "SMESH_MesherHelper.hxx"
32 #include "SMESH_subMesh.hxx"
33 #include "SMESH_subMeshEventListener.hxx"
34 #include "StdMeshers_CartesianParameters3D.hxx"
35
36 #include <utilities.h>
37 #include <Utils_ExceptHandlers.hxx>
38 #include <Basics_OCCTVersion.hxx>
39
40 #include <GEOMUtils.hxx>
41
42 #include <BRepAdaptor_Curve.hxx>
43 #include <BRepAdaptor_Surface.hxx>
44 #include <BRepBndLib.hxx>
45 #include <BRepBuilderAPI_Copy.hxx>
46 #include <BRepBuilderAPI_MakeFace.hxx>
47 #include <BRepTools.hxx>
48 #include <BRep_Builder.hxx>
49 #include <BRep_Tool.hxx>
50 #include <Bnd_B3d.hxx>
51 #include <Bnd_Box.hxx>
52 #include <ElSLib.hxx>
53 #include <GCPnts_UniformDeflection.hxx>
54 #include <Geom2d_BSplineCurve.hxx>
55 #include <Geom2d_BezierCurve.hxx>
56 #include <Geom2d_TrimmedCurve.hxx>
57 #include <GeomAPI_ProjectPointOnSurf.hxx>
58 #include <GeomLib.hxx>
59 #include <Geom_BSplineCurve.hxx>
60 #include <Geom_BSplineSurface.hxx>
61 #include <Geom_BezierCurve.hxx>
62 #include <Geom_BezierSurface.hxx>
63 #include <Geom_RectangularTrimmedSurface.hxx>
64 #include <Geom_TrimmedCurve.hxx>
65 #include <IntAna_IntConicQuad.hxx>
66 #include <IntAna_IntLinTorus.hxx>
67 #include <IntAna_Quadric.hxx>
68 #include <IntCurveSurface_TransitionOnCurve.hxx>
69 #include <IntCurvesFace_Intersector.hxx>
70 #include <Poly_Triangulation.hxx>
71 #include <Precision.hxx>
72 #include <TopExp.hxx>
73 #include <TopExp_Explorer.hxx>
74 #include <TopLoc_Location.hxx>
75 #include <TopTools_MapOfShape.hxx>
76 #include <TopoDS.hxx>
77 #include <TopoDS_Compound.hxx>
78 #include <TopoDS_Face.hxx>
79 #include <TopoDS_TShape.hxx>
80 #include <gp_Cone.hxx>
81 #include <gp_Cylinder.hxx>
82 #include <gp_Lin.hxx>
83 #include <gp_Pln.hxx>
84 #include <gp_Pnt2d.hxx>
85 #include <gp_Sphere.hxx>
86 #include <gp_Torus.hxx>
87
88 #include <limits>
89
90 //#undef WITH_TBB
91 #ifdef WITH_TBB
92 #include <tbb/parallel_for.h>
93 //#include <tbb/enumerable_thread_specific.h>
94 #endif
95
96 using namespace std;
97
98 #ifdef _DEBUG_
99 //#define _MY_DEBUG_
100 #endif
101
102 #if OCC_VERSION_LARGE <= 0x06050300
103 // workaround is required only for OCCT6.5.3 and older (see OCC22809)
104 #define ELLIPSOLID_WORKAROUND
105 #endif
106
107 #ifdef ELLIPSOLID_WORKAROUND
108 #include <BRepIntCurveSurface_Inter.hxx>
109 #include <BRepTopAdaptor_TopolTool.hxx>
110 #include <BRepAdaptor_HSurface.hxx>
111 #endif
112
113 //=============================================================================
114 /*!
115  * Constructor
116  */
117 //=============================================================================
118
119 StdMeshers_Cartesian_3D::StdMeshers_Cartesian_3D(int hypId, int studyId, SMESH_Gen * gen)
120   :SMESH_3D_Algo(hypId, studyId, gen)
121 {
122   _name = "Cartesian_3D";
123   _shapeType = (1 << TopAbs_SOLID);       // 1 bit /shape type
124   _compatibleHypothesis.push_back("CartesianParameters3D");
125
126   _onlyUnaryInput = false;          // to mesh all SOLIDs at once
127   _requireDiscreteBoundary = false; // 2D mesh not needed
128   _supportSubmeshes = false;        // do not use any existing mesh
129 }
130
131 //=============================================================================
132 /*!
133  * Check presence of a hypothesis
134  */
135 //=============================================================================
136
137 bool StdMeshers_Cartesian_3D::CheckHypothesis (SMESH_Mesh&          aMesh,
138                                                const TopoDS_Shape&  aShape,
139                                                Hypothesis_Status&   aStatus)
140 {
141   aStatus = SMESH_Hypothesis::HYP_MISSING;
142
143   const list<const SMESHDS_Hypothesis*>& hyps = GetUsedHypothesis(aMesh, aShape);
144   list <const SMESHDS_Hypothesis* >::const_iterator h = hyps.begin();
145   if ( h == hyps.end())
146   {
147     return false;
148   }
149
150   for ( ; h != hyps.end(); ++h )
151   {
152     if (( _hyp = dynamic_cast<const StdMeshers_CartesianParameters3D*>( *h )))
153     {
154       aStatus = _hyp->IsDefined() ? HYP_OK : HYP_BAD_PARAMETER;
155       break;
156     }
157   }
158
159   return aStatus == HYP_OK;
160 }
161
162 namespace
163 {
164   typedef int TGeomID;
165
166   //=============================================================================
167   // Definitions of internal utils
168   // --------------------------------------------------------------------------
169   enum Transition {
170     Trans_TANGENT = IntCurveSurface_Tangent,
171     Trans_IN      = IntCurveSurface_In,
172     Trans_OUT     = IntCurveSurface_Out,
173     Trans_APEX
174   };
175   // --------------------------------------------------------------------------
176   /*!
177    * \brief Common data of any intersection between a Grid and a shape
178    */
179   struct B_IntersectPoint
180   {
181     mutable const SMDS_MeshNode* _node;
182     mutable vector< TGeomID >    _faceIDs;
183
184     B_IntersectPoint(): _node(NULL) {}
185     void Add( const vector< TGeomID >& fIDs, const SMDS_MeshNode* n=0 ) const;
186     int HasCommonFace( const B_IntersectPoint * other, int avoidFace=-1 ) const;
187     bool IsOnFace( int faceID ) const;
188     virtual ~B_IntersectPoint() {}
189   };
190   // --------------------------------------------------------------------------
191   /*!
192    * \brief Data of intersection between a GridLine and a TopoDS_Face
193    */
194   struct F_IntersectPoint : public B_IntersectPoint
195   {
196     double             _paramOnLine;
197     mutable Transition _transition;
198     mutable size_t     _indexOnLine;
199
200     bool operator< ( const F_IntersectPoint& o ) const { return _paramOnLine < o._paramOnLine; }
201   };
202   // --------------------------------------------------------------------------
203   /*!
204    * \brief Data of intersection between GridPlanes and a TopoDS_EDGE
205    */
206   struct E_IntersectPoint : public B_IntersectPoint
207   {
208     gp_Pnt  _point;
209     double  _uvw[3];
210     TGeomID _shapeID;
211   };
212   // --------------------------------------------------------------------------
213   /*!
214    * \brief A line of the grid and its intersections with 2D geometry
215    */
216   struct GridLine
217   {
218     gp_Lin _line;
219     double _length; // line length
220     multiset< F_IntersectPoint > _intPoints;
221
222     void RemoveExcessIntPoints( const double tol );
223     bool GetIsOutBefore( multiset< F_IntersectPoint >::iterator ip, bool prevIsOut );
224   };
225   // --------------------------------------------------------------------------
226   /*!
227    * \brief Planes of the grid used to find intersections of an EDGE with a hexahedron
228    */
229   struct GridPlanes
230   {
231     gp_XYZ           _zNorm;
232     vector< gp_XYZ > _origins; // origin points of all planes in one direction
233     vector< double > _zProjs;  // projections of origins to _zNorm
234   };
235   // --------------------------------------------------------------------------
236   /*!
237    * \brief Iterator on the parallel grid lines of one direction
238    */
239   struct LineIndexer
240   {
241     size_t _size  [3];
242     size_t _curInd[3];
243     size_t _iVar1, _iVar2, _iConst;
244     string _name1, _name2, _nameConst;
245     LineIndexer() {}
246     LineIndexer( size_t sz1, size_t sz2, size_t sz3,
247                  size_t iv1, size_t iv2, size_t iConst,
248                  const string& nv1, const string& nv2, const string& nConst )
249     {
250       _size[0] = sz1; _size[1] = sz2; _size[2] = sz3;
251       _curInd[0] = _curInd[1] = _curInd[2] = 0;
252       _iVar1 = iv1; _iVar2 = iv2; _iConst = iConst; 
253       _name1 = nv1; _name2 = nv2; _nameConst = nConst;
254     }
255
256     size_t I() const { return _curInd[0]; }
257     size_t J() const { return _curInd[1]; }
258     size_t K() const { return _curInd[2]; }
259     void SetIJK( size_t i, size_t j, size_t k )
260     {
261       _curInd[0] = i; _curInd[1] = j; _curInd[2] = k;
262     }
263     void operator++()
264     {
265       if ( ++_curInd[_iVar1] == _size[_iVar1] )
266         _curInd[_iVar1] = 0, ++_curInd[_iVar2];
267     }
268     bool More() const { return _curInd[_iVar2] < _size[_iVar2]; }
269     size_t LineIndex   () const { return _curInd[_iVar1] + _curInd[_iVar2]* _size[_iVar1]; }
270     size_t LineIndex10 () const { return (_curInd[_iVar1] + 1 ) + _curInd[_iVar2]* _size[_iVar1]; }
271     size_t LineIndex01 () const { return _curInd[_iVar1] + (_curInd[_iVar2] + 1 )* _size[_iVar1]; }
272     size_t LineIndex11 () const { return (_curInd[_iVar1] + 1 ) + (_curInd[_iVar2] + 1 )* _size[_iVar1]; }
273     void SetIndexOnLine (size_t i)  { _curInd[ _iConst ] = i; }
274     size_t NbLines() const { return _size[_iVar1] * _size[_iVar2]; }
275   };
276   // --------------------------------------------------------------------------
277   /*!
278    * \brief Container of GridLine's
279    */
280   struct Grid
281   {
282     vector< double >   _coords[3]; // coordinates of grid nodes
283     gp_XYZ             _axes  [3]; // axis directions
284     vector< GridLine > _lines [3]; //    in 3 directions
285     double             _tol, _minCellSize;
286     gp_XYZ             _origin;
287     gp_Mat             _invB; // inverted basis of _axes
288
289     vector< const SMDS_MeshNode* >    _nodes; // mesh nodes at grid nodes
290     vector< const F_IntersectPoint* > _gridIntP; // grid node intersection with geometry
291
292     list< E_IntersectPoint >          _edgeIntP; // intersections with EDGEs
293     TopTools_IndexedMapOfShape        _shapes;
294
295     SMESH_MesherHelper*               _helper;
296
297     size_t CellIndex( size_t i, size_t j, size_t k ) const
298     {
299       return i + j*(_coords[0].size()-1) + k*(_coords[0].size()-1)*(_coords[1].size()-1);
300     }
301     size_t NodeIndex( size_t i, size_t j, size_t k ) const
302     {
303       return i + j*_coords[0].size() + k*_coords[0].size()*_coords[1].size();
304     }
305     size_t NodeIndexDX() const { return 1; }
306     size_t NodeIndexDY() const { return _coords[0].size(); }
307     size_t NodeIndexDZ() const { return _coords[0].size() * _coords[1].size(); }
308
309     LineIndexer GetLineIndexer(size_t iDir) const;
310
311     void SetCoordinates(const vector<double>& xCoords,
312                         const vector<double>& yCoords,
313                         const vector<double>& zCoords,
314                         const double*         axesDirs,
315                         const Bnd_Box&        bndBox );
316     void ComputeUVW(const gp_XYZ& p, double uvw[3]);
317     void ComputeNodes(SMESH_MesherHelper& helper);
318   };
319 #ifdef ELLIPSOLID_WORKAROUND
320   // --------------------------------------------------------------------------
321   /*!
322    * \brief struct temporary replacing IntCurvesFace_Intersector until
323    *        OCCT bug 0022809 is fixed
324    *        http://tracker.dev.opencascade.org/view.php?id=22809
325    */
326   struct TMP_IntCurvesFace_Intersector
327   {
328     BRepAdaptor_Surface                       _surf;
329     double                                    _tol;
330     BRepIntCurveSurface_Inter                 _intcs;
331     vector<IntCurveSurface_IntersectionPoint> _points;
332     BRepTopAdaptor_TopolTool                  _clsf;
333
334     TMP_IntCurvesFace_Intersector(const TopoDS_Face& face, const double tol)
335       :_surf( face ), _tol( tol ), _clsf( new BRepAdaptor_HSurface(_surf) ) {}
336     Bnd_Box Bounding() const { Bnd_Box b; BRepBndLib::Add (_surf.Face(), b); return b; }
337     void Perform( const gp_Lin& line, const double w0, const double w1 )
338     {
339       _points.clear();
340       for ( _intcs.Init( _surf.Face(), line, _tol ); _intcs.More(); _intcs.Next() )
341         if ( w0 <= _intcs.W() && _intcs.W() <= w1 )
342           _points.push_back( _intcs.Point() );
343     }
344     bool IsDone() const { return true; }
345     int  NbPnt()  const { return _points.size(); }
346     IntCurveSurface_TransitionOnCurve Transition( const int i ) const { return _points[ i-1 ].Transition(); }
347     double       WParameter( const int i ) const { return _points[ i-1 ].W(); }
348     TopAbs_State ClassifyUVPoint(const gp_Pnt2d& p) { return _clsf.Classify( p, _tol ); }
349   };
350 #define __IntCurvesFace_Intersector TMP_IntCurvesFace_Intersector
351 #else
352 #define __IntCurvesFace_Intersector IntCurvesFace_Intersector
353 #endif
354   // --------------------------------------------------------------------------
355   /*!
356    * \brief Intersector of TopoDS_Face with all GridLine's
357    */
358   struct FaceGridIntersector
359   {
360     TopoDS_Face _face;
361     TGeomID     _faceID;
362     Grid*       _grid;
363     Bnd_Box     _bndBox;
364     __IntCurvesFace_Intersector* _surfaceInt;
365     vector< std::pair< GridLine*, F_IntersectPoint > > _intersections;
366
367     FaceGridIntersector(): _grid(0), _surfaceInt(0) {}
368     void Intersect();
369
370     void StoreIntersections()
371     {
372       for ( size_t i = 0; i < _intersections.size(); ++i )
373       {
374         multiset< F_IntersectPoint >::iterator ip = 
375           _intersections[i].first->_intPoints.insert( _intersections[i].second );
376         ip->_faceIDs.reserve( 1 );
377         ip->_faceIDs.push_back( _faceID );
378       }
379     }
380     const Bnd_Box& GetFaceBndBox()
381     {
382       GetCurveFaceIntersector();
383       return _bndBox;
384     }
385     __IntCurvesFace_Intersector* GetCurveFaceIntersector()
386     {
387       if ( !_surfaceInt )
388       {
389         _surfaceInt = new __IntCurvesFace_Intersector( _face, Precision::PConfusion() );
390         _bndBox     = _surfaceInt->Bounding();
391         if ( _bndBox.IsVoid() )
392           BRepBndLib::Add (_face, _bndBox);
393       }
394       return _surfaceInt;
395     }
396     bool IsThreadSafe(set< const Standard_Transient* >& noSafeTShapes) const;
397   };
398   // --------------------------------------------------------------------------
399   /*!
400    * \brief Intersector of a surface with a GridLine
401    */
402   struct FaceLineIntersector
403   {
404     double      _tol;
405     double      _u, _v, _w; // params on the face and the line
406     Transition  _transition; // transition of at intersection (see IntCurveSurface.cdl)
407     Transition  _transIn, _transOut; // IN and OUT transitions depending of face orientation
408
409     gp_Pln      _plane;
410     gp_Cylinder _cylinder;
411     gp_Cone     _cone;
412     gp_Sphere   _sphere;
413     gp_Torus    _torus;
414     __IntCurvesFace_Intersector* _surfaceInt;
415
416     vector< F_IntersectPoint > _intPoints;
417
418     void IntersectWithPlane   (const GridLine& gridLine);
419     void IntersectWithCylinder(const GridLine& gridLine);
420     void IntersectWithCone    (const GridLine& gridLine);
421     void IntersectWithSphere  (const GridLine& gridLine);
422     void IntersectWithTorus   (const GridLine& gridLine);
423     void IntersectWithSurface (const GridLine& gridLine);
424
425     bool UVIsOnFace() const;
426     void addIntPoint(const bool toClassify=true);
427     bool isParamOnLineOK( const double linLength )
428     {
429       return -_tol < _w && _w < linLength + _tol;
430     }
431     FaceLineIntersector():_surfaceInt(0) {}
432     ~FaceLineIntersector() { if (_surfaceInt ) delete _surfaceInt; _surfaceInt = 0; }
433   };
434   // --------------------------------------------------------------------------
435   /*!
436    * \brief Class representing topology of the hexahedron and creating a mesh
437    *        volume basing on analysis of hexahedron intersection with geometry
438    */
439   class Hexahedron
440   {
441     // --------------------------------------------------------------------------------
442     struct _Face;
443     struct _Link;
444     // --------------------------------------------------------------------------------
445     struct _Node //!< node either at a hexahedron corner or at intersection
446     {
447       const SMDS_MeshNode*    _node; // mesh node at hexahedron corner
448       const B_IntersectPoint* _intPoint;
449       const _Face*            _usedInFace;
450
451       _Node(const SMDS_MeshNode* n=0, const B_IntersectPoint* ip=0)
452         :_node(n), _intPoint(ip), _usedInFace(0) {} 
453       const SMDS_MeshNode*    Node() const
454       { return ( _intPoint && _intPoint->_node ) ? _intPoint->_node : _node; }
455       //const F_IntersectPoint* FaceIntPnt() const
456       //{ return static_cast< const F_IntersectPoint* >( _intPoint ); }
457       const E_IntersectPoint* EdgeIntPnt() const
458       { return static_cast< const E_IntersectPoint* >( _intPoint ); }
459       bool IsUsedInFace( const _Face* polygon = 0 )
460       {
461         return polygon ? ( _usedInFace == polygon ) : bool( _usedInFace );
462       }
463       void Add( const E_IntersectPoint* ip )
464       {
465         if ( !_intPoint ) {
466           _intPoint = ip;
467         }
468         else if ( !_intPoint->_node ) {
469           ip->Add( _intPoint->_faceIDs );
470           _intPoint = ip;
471         }
472         else  {
473           _intPoint->Add( ip->_faceIDs );
474         }
475       }
476       int IsLinked( const B_IntersectPoint* other,
477                     int                     avoidFace=-1 ) const // returns id of a common face
478       {
479         return _intPoint ? _intPoint->HasCommonFace( other, avoidFace ) : 0;
480       }
481       bool IsOnFace( int faceID ) const // returns true if faceID is found
482       {
483         return _intPoint ? _intPoint->IsOnFace( faceID ) : false;
484       }
485       gp_Pnt Point() const
486       {
487         if ( const SMDS_MeshNode* n = Node() )
488           return SMESH_TNodeXYZ( n );
489         if ( const E_IntersectPoint* eip =
490              dynamic_cast< const E_IntersectPoint* >( _intPoint ))
491           return eip->_point;
492         return gp_Pnt( 1e100, 0, 0 );
493       }
494     };
495     // --------------------------------------------------------------------------------
496     struct _Link // link connecting two _Node's
497     {
498       _Node* _nodes[2];
499       _Face* _faces[2]; // polygons sharing a link
500       vector< const F_IntersectPoint* > _fIntPoints; // GridLine intersections with FACEs
501       vector< _Node* >                  _fIntNodes;   // _Node's at _fIntPoints
502       vector< _Link >                   _splits;
503       _Link() { _faces[0] = 0; }
504     };
505     // --------------------------------------------------------------------------------
506     struct _OrientedLink
507     {
508       _Link* _link;
509       bool   _reverse;
510       _OrientedLink( _Link* link=0, bool reverse=false ): _link(link), _reverse(reverse) {}
511       void Reverse() { _reverse = !_reverse; }
512       int NbResultLinks() const { return _link->_splits.size(); }
513       _OrientedLink ResultLink(int i) const
514       {
515         return _OrientedLink(&_link->_splits[_reverse ? NbResultLinks()-i-1 : i],_reverse);
516       }
517       _Node* FirstNode() const { return _link->_nodes[ _reverse ]; }
518       _Node* LastNode()  const { return _link->_nodes[ !_reverse ]; }
519       operator bool() const { return _link; }
520       vector< TGeomID > GetNotUsedFace(const set<TGeomID>& usedIDs ) const // returns supporting FACEs
521       {
522         vector< TGeomID > faces;
523         const B_IntersectPoint *ip0, *ip1;
524         if (( ip0 = _link->_nodes[0]->_intPoint ) &&
525             ( ip1 = _link->_nodes[1]->_intPoint ))
526         {
527           for ( size_t i = 0; i < ip0->_faceIDs.size(); ++i )
528             if ( ip1->IsOnFace ( ip0->_faceIDs[i] ) &&
529                  !usedIDs.count( ip0->_faceIDs[i] ) )
530               faces.push_back( ip0->_faceIDs[i] );
531         }
532         return faces;
533       }
534       bool HasEdgeNodes() const
535       {
536         return ( dynamic_cast< const E_IntersectPoint* >( _link->_nodes[0]->_intPoint ) ||
537                  dynamic_cast< const E_IntersectPoint* >( _link->_nodes[1]->_intPoint ));
538       }
539       int NbFaces() const
540       {
541         return !_link->_faces[0] ? 0 : 1 + bool( _link->_faces[1] );
542       }
543       void AddFace( _Face* f )
544       {
545         if ( _link->_faces[0] )
546         {
547           _link->_faces[1] = f;
548         }
549         else
550         {
551           _link->_faces[0] = f;
552           _link->_faces[1] = 0;
553         }
554       }
555       void RemoveFace( _Face* f )
556       {
557         if ( !_link->_faces[0] ) return;
558
559         if ( _link->_faces[1] == f )
560         {
561           _link->_faces[1] = 0;
562         }
563         else if ( _link->_faces[0] == f )
564         {
565           _link->_faces[0];
566           if ( _link->_faces[1] )
567           {
568             _link->_faces[0] = _link->_faces[1];
569             _link->_faces[1] = 0;
570           }
571         }
572       }
573     };
574     // --------------------------------------------------------------------------------
575     struct _Face
576     {
577       vector< _OrientedLink > _links;       // links on GridLine's
578       vector< _Link >         _polyLinks;   // links added to close a polygonal face
579       vector< _Node* >        _eIntNodes;   // nodes at intersection with EDGEs
580       bool isPolyLink( const _OrientedLink& ol )
581       {
582         return _polyLinks.empty() ? false :
583           ( &_polyLinks[0] <= ol._link &&  ol._link <= &_polyLinks.back() );
584       }
585     };
586     // --------------------------------------------------------------------------------
587     struct _volumeDef // holder of nodes of a volume mesh element
588     {
589       vector< _Node* > _nodes;
590       vector< int >    _quantities;
591       typedef boost::shared_ptr<_volumeDef> Ptr;
592       void set( const vector< _Node* >& nodes,
593                 const vector< int >&    quant = vector< int >() )
594       { _nodes = nodes; _quantities = quant; }
595     };
596
597     // topology of a hexahedron
598     int   _nodeShift[8];
599     _Node _hexNodes [8];
600     _Link _hexLinks [12];
601     _Face _hexQuads [6];
602
603     // faces resulted from hexahedron intersection
604     vector< _Face > _polygons;
605
606     // intresections with EDGEs
607     vector< const E_IntersectPoint* > _eIntPoints;
608
609     // additional nodes created at intersection points
610     vector< _Node > _intNodes;
611
612     // nodes inside the hexahedron (at VERTEXes)
613     vector< _Node* > _vIntNodes;
614
615     // computed volume elements
616     //vector< _volumeDef::Ptr > _volumeDefs;
617     _volumeDef _volumeDefs;
618
619     Grid*       _grid;
620     double      _sizeThreshold, _sideLength[3];
621     int         _nbCornerNodes, _nbFaceIntNodes, _nbBndNodes;
622     int         _origNodeInd; // index of _hexNodes[0] node within the _grid
623     size_t      _i,_j,_k;
624
625   public:
626     Hexahedron(const double sizeThreshold, Grid* grid);
627     int MakeElements(SMESH_MesherHelper&                      helper,
628                      const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap);
629     void ComputeElements();
630     void Init() { init( _i, _j, _k ); }
631
632   private:
633     Hexahedron(const Hexahedron& other );
634     void init( size_t i, size_t j, size_t k );
635     void init( size_t i );
636     void addEdges(SMESH_MesherHelper&                      helper,
637                   vector< Hexahedron* >&                   intersectedHex,
638                   const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap);
639     gp_Pnt findIntPoint( double u1, double proj1, double u2, double proj2,
640                          double proj, BRepAdaptor_Curve& curve,
641                          const gp_XYZ& axis, const gp_XYZ& origin );
642     int  getEntity( const E_IntersectPoint* ip, int* facets, int& sub );
643     bool addIntersection( const E_IntersectPoint& ip,
644                           vector< Hexahedron* >&  hexes,
645                           int ijk[], int dIJK[] );
646     bool findChain( _Node* n1, _Node* n2, _Face& quad, vector<_Node*>& chainNodes );
647     bool closePolygon( _Face* polygon, vector<_Node*>& chainNodes ) const;
648     int  addElements(SMESH_MesherHelper& helper);
649     bool is1stNodeOut( _Link& link ) const;
650     bool isInHole() const;
651     bool checkPolyhedronSize() const;
652     bool addHexa ();
653     bool addTetra();
654     bool addPenta();
655     bool addPyra ();
656     bool debugDumpLink( _Link* link );
657     _Node* FindEqualNode( vector< _Node* >&       nodes,
658                           const E_IntersectPoint* ip,
659                           const double            tol2 )
660     {
661       for ( size_t i = 0; i < nodes.size(); ++i )
662         if ( nodes[i]->EdgeIntPnt() == ip ||
663              nodes[i]->Point().SquareDistance( ip->_point ) <= tol2 )
664           return nodes[i];
665       return 0;
666     }
667   };
668
669 #ifdef WITH_TBB
670   // --------------------------------------------------------------------------
671   /*!
672    * \brief Hexahedron computing volumes in one thread
673    */
674   struct ParallelHexahedron
675   {
676     vector< Hexahedron* >& _hexVec;
677     vector<int>&           _index;
678     ParallelHexahedron( vector< Hexahedron* >& hv, vector<int>& ind): _hexVec(hv), _index(ind) {}
679     void operator() ( const tbb::blocked_range<size_t>& r ) const
680     {
681       for ( size_t i = r.begin(); i != r.end(); ++i )
682         if ( Hexahedron* hex = _hexVec[ _index[i]] )
683           hex->ComputeElements();
684     }
685   };
686   // --------------------------------------------------------------------------
687   /*!
688    * \brief Structure intersecting certain nb of faces with GridLine's in one thread
689    */
690   struct ParallelIntersector
691   {
692     vector< FaceGridIntersector >& _faceVec;
693     ParallelIntersector( vector< FaceGridIntersector >& faceVec): _faceVec(faceVec){}
694     void operator() ( const tbb::blocked_range<size_t>& r ) const
695     {
696       for ( size_t i = r.begin(); i != r.end(); ++i )
697         _faceVec[i].Intersect();
698     }
699   };
700 #endif
701
702   //=============================================================================
703   // Implementation of internal utils
704   //=============================================================================
705   /*!
706    * \brief adjust \a i to have \a val between values[i] and values[i+1]
707    */
708   inline void locateValue( int & i, double val, const vector<double>& values,
709                            int& di, double tol )
710   {
711     //val += values[0]; // input \a val is measured from 0.
712     if ( i > values.size()-2 )
713       i = values.size()-2;
714     else
715       while ( i+2 < values.size() && val > values[ i+1 ])
716         ++i;
717     while ( i > 0 && val < values[ i ])
718       --i;
719
720     if ( i > 0 && val - values[ i ] < tol )
721       di = -1;
722     else if ( i+2 < values.size() && values[ i+1 ] - val < tol )
723       di = 1;
724     else
725       di = 0;
726   }
727   //=============================================================================
728   /*
729    * Remove coincident intersection points
730    */
731   void GridLine::RemoveExcessIntPoints( const double tol )
732   {
733     if ( _intPoints.size() < 2 ) return;
734
735     set< Transition > tranSet;
736     multiset< F_IntersectPoint >::iterator ip1, ip2 = _intPoints.begin();
737     while ( ip2 != _intPoints.end() )
738     {
739       tranSet.clear();
740       ip1 = ip2++;
741       while ( ip2 != _intPoints.end() && ip2->_paramOnLine - ip1->_paramOnLine <= tol )
742       {
743         tranSet.insert( ip1->_transition );
744         tranSet.insert( ip2->_transition );
745         ip2->Add( ip1->_faceIDs );
746         _intPoints.erase( ip1 );
747         ip1 = ip2++;
748       }
749       if ( tranSet.size() > 1 ) // points with different transition coincide
750       {
751         bool isIN  = tranSet.count( Trans_IN );
752         bool isOUT = tranSet.count( Trans_OUT );
753         if ( isIN && isOUT )
754           (*ip1)._transition = Trans_TANGENT;
755         else
756           (*ip1)._transition = isIN ? Trans_IN : Trans_OUT;
757       }
758     }
759   }
760   //================================================================================
761   /*
762    * Return "is OUT" state for nodes before the given intersection point
763    */
764   bool GridLine::GetIsOutBefore( multiset< F_IntersectPoint >::iterator ip, bool prevIsOut )
765   {
766     if ( ip->_transition == Trans_IN )
767       return true;
768     if ( ip->_transition == Trans_OUT )
769       return false;
770     if ( ip->_transition == Trans_APEX )
771     {
772       // singularity point (apex of a cone)
773       if ( _intPoints.size() == 1 || ip == _intPoints.begin() )
774         return true;
775       multiset< F_IntersectPoint >::iterator ipBef = ip, ipAft = ++ip;
776       if ( ipAft == _intPoints.end() )
777         return false;
778       --ipBef;
779       if ( ipBef->_transition != ipAft->_transition )
780         return ( ipBef->_transition == Trans_OUT );
781       return ( ipBef->_transition != Trans_OUT );
782     }
783     // _transition == Trans_TANGENT
784     return !prevIsOut;
785   }
786   //================================================================================
787   /*
788    * Adds face IDs
789    */
790   void B_IntersectPoint::Add( const vector< TGeomID >& fIDs,
791                               const SMDS_MeshNode*     n) const
792   {
793     if ( _faceIDs.empty() )
794       _faceIDs = fIDs;
795     else
796       for ( size_t i = 0; i < fIDs.size(); ++i )
797       {
798         vector< TGeomID >::iterator it =
799           std::find( _faceIDs.begin(), _faceIDs.end(), fIDs[i] );
800         if ( it == _faceIDs.end() )
801           _faceIDs.push_back( fIDs[i] );
802       }
803     if ( !_node )
804       _node = n;
805   }
806   //================================================================================
807   /*
808    * Returns index of a common face if any, else zero
809    */
810   int B_IntersectPoint::HasCommonFace( const B_IntersectPoint * other, int avoidFace ) const
811   {
812     if ( other )
813       for ( size_t i = 0; i < other->_faceIDs.size(); ++i )
814         if ( avoidFace != other->_faceIDs[i] &&
815              IsOnFace   ( other->_faceIDs[i] ))
816           return other->_faceIDs[i];
817     return 0;
818   }
819   //================================================================================
820   /*
821    * Returns \c true if \a faceID in in this->_faceIDs
822    */
823   bool B_IntersectPoint::IsOnFace( int faceID ) const // returns true if faceID is found
824   {
825     vector< TGeomID >::const_iterator it =
826       std::find( _faceIDs.begin(), _faceIDs.end(), faceID );
827     return ( it != _faceIDs.end() );
828   }
829   //================================================================================
830   /*
831    * Return an iterator on GridLine's in a given direction
832    */
833   LineIndexer Grid::GetLineIndexer(size_t iDir) const
834   {
835     const size_t indices[] = { 1,2,0, 0,2,1, 0,1,2 };
836     const string s      [] = { "X", "Y", "Z" };
837     LineIndexer li( _coords[0].size(),  _coords[1].size(),    _coords[2].size(),
838                     indices[iDir*3],    indices[iDir*3+1],    indices[iDir*3+2],
839                     s[indices[iDir*3]], s[indices[iDir*3+1]], s[indices[iDir*3+2]]);
840     return li;
841   }
842   //=============================================================================
843   /*
844    * Creates GridLine's of the grid
845    */
846   void Grid::SetCoordinates(const vector<double>& xCoords,
847                             const vector<double>& yCoords,
848                             const vector<double>& zCoords,
849                             const double*         axesDirs,
850                             const Bnd_Box&        shapeBox)
851   {
852     _coords[0] = xCoords;
853     _coords[1] = yCoords;
854     _coords[2] = zCoords;
855
856     _axes[0].SetCoord( axesDirs[0],
857                        axesDirs[1],
858                        axesDirs[2]);
859     _axes[1].SetCoord( axesDirs[3],
860                        axesDirs[4],
861                        axesDirs[5]);
862     _axes[2].SetCoord( axesDirs[6],
863                        axesDirs[7],
864                        axesDirs[8]);
865     _axes[0].Normalize();
866     _axes[1].Normalize();
867     _axes[2].Normalize();
868
869     _invB.SetCols( _axes[0], _axes[1], _axes[2] );
870     _invB.Invert();
871
872     // compute tolerance
873     _minCellSize = Precision::Infinite();
874     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
875     {
876       for ( size_t i = 1; i < _coords[ iDir ].size(); ++i )
877       {
878         double cellLen = _coords[ iDir ][ i ] - _coords[ iDir ][ i-1 ];
879         if ( cellLen < _minCellSize )
880           _minCellSize = cellLen;
881       }
882     }
883     if ( _minCellSize < Precision::Confusion() )
884       throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
885                                 SMESH_Comment("Too small cell size: ") << _minCellSize );
886     _tol = _minCellSize / 1000.;
887
888     // attune grid extremities to shape bounding box
889
890     double sP[6]; // aXmin, aYmin, aZmin, aXmax, aYmax, aZmax
891     shapeBox.Get(sP[0],sP[1],sP[2],sP[3],sP[4],sP[5]);
892     double* cP[6] = { &_coords[0].front(), &_coords[1].front(), &_coords[2].front(),
893                       &_coords[0].back(),  &_coords[1].back(),  &_coords[2].back() };
894     for ( int i = 0; i < 6; ++i )
895       if ( fabs( sP[i] - *cP[i] ) < _tol )
896         *cP[i] = sP[i];// + _tol/1000. * ( i < 3 ? +1 : -1 );
897
898     for ( int iDir = 0; iDir < 3; ++iDir )
899     {
900       if ( _coords[iDir][0] - sP[iDir] > _tol )
901       {
902         _minCellSize = Min( _minCellSize, _coords[iDir][0] - sP[iDir] );
903         _coords[iDir].insert( _coords[iDir].begin(), sP[iDir] + _tol/1000.);
904       }
905       if ( sP[iDir+3] - _coords[iDir].back() > _tol  )
906       {
907         _minCellSize = Min( _minCellSize, sP[iDir+3] - _coords[iDir].back() );
908         _coords[iDir].push_back( sP[iDir+3] - _tol/1000.);
909       }
910     }
911     _tol = _minCellSize / 1000.;
912
913     _origin = ( _coords[0][0] * _axes[0] +
914                 _coords[1][0] * _axes[1] +
915                 _coords[2][0] * _axes[2] );
916
917     // create lines
918     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
919     {
920       LineIndexer li = GetLineIndexer( iDir );
921       _lines[iDir].resize( li.NbLines() );
922       double len = _coords[ iDir ].back() - _coords[iDir].front();
923       for ( ; li.More(); ++li )
924       {
925         GridLine& gl = _lines[iDir][ li.LineIndex() ];
926         gl._line.SetLocation( _coords[0][li.I()] * _axes[0] +
927                               _coords[1][li.J()] * _axes[1] +
928                               _coords[2][li.K()] * _axes[2] );
929         gl._line.SetDirection( _axes[ iDir ]);
930         gl._length = len;
931       }
932     }
933   }
934   //================================================================================
935   /*
936    * Computes coordinates of a point in the grid CS
937    */
938   void Grid::ComputeUVW(const gp_XYZ& P, double UVW[3])
939   {
940     gp_XYZ p = P * _invB;
941     p.Coord( UVW[0], UVW[1], UVW[2] );
942   }
943   //================================================================================
944   /*
945    * Creates all nodes
946    */
947   void Grid::ComputeNodes(SMESH_MesherHelper& helper)
948   {
949     // state of each node of the grid relative to the geometry
950     const size_t nbGridNodes = _coords[0].size() * _coords[1].size() * _coords[2].size();
951     vector< bool > isNodeOut( nbGridNodes, false );
952     _nodes.resize( nbGridNodes, 0 );
953     _gridIntP.resize( nbGridNodes, NULL );
954
955     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
956     {
957       LineIndexer li = GetLineIndexer( iDir );
958
959       // find out a shift of node index while walking along a GridLine in this direction
960       li.SetIndexOnLine( 0 );
961       size_t nIndex0 = NodeIndex( li.I(), li.J(), li.K() );
962       li.SetIndexOnLine( 1 );
963       const size_t nShift = NodeIndex( li.I(), li.J(), li.K() ) - nIndex0;
964       
965       const vector<double> & coords = _coords[ iDir ];
966       for ( ; li.More(); ++li ) // loop on lines in iDir
967       {
968         li.SetIndexOnLine( 0 );
969         nIndex0 = NodeIndex( li.I(), li.J(), li.K() );
970
971         GridLine& line = _lines[ iDir ][ li.LineIndex() ];
972         const gp_XYZ lineLoc = line._line.Location().XYZ();
973         const gp_XYZ lineDir = line._line.Direction().XYZ();
974         line.RemoveExcessIntPoints( _tol );
975         multiset< F_IntersectPoint >& intPnts = line._intPoints;
976         multiset< F_IntersectPoint >::iterator ip = intPnts.begin();
977
978         bool isOut = true;
979         const double* nodeCoord = & coords[0];
980         const double* coord0    = nodeCoord;
981         const double* coordEnd  = coord0 + coords.size();
982         double nodeParam = 0;
983         for ( ; ip != intPnts.end(); ++ip )
984         {
985           // set OUT state or just skip IN nodes before ip
986           if ( nodeParam < ip->_paramOnLine - _tol )
987           {
988             isOut = line.GetIsOutBefore( ip, isOut );
989
990             while ( nodeParam < ip->_paramOnLine - _tol )
991             {
992               if ( isOut )
993                 isNodeOut[ nIndex0 + nShift * ( nodeCoord-coord0 ) ] = isOut;
994               if ( ++nodeCoord <  coordEnd )
995                 nodeParam = *nodeCoord - *coord0;
996               else
997                 break;
998             }
999             if ( nodeCoord == coordEnd ) break;
1000           }
1001           // create a mesh node on a GridLine at ip if it does not coincide with a grid node
1002           if ( nodeParam > ip->_paramOnLine + _tol )
1003           {
1004             // li.SetIndexOnLine( 0 );
1005             // double xyz[3] = { _coords[0][ li.I() ], _coords[1][ li.J() ], _coords[2][ li.K() ]};
1006             // xyz[ li._iConst ] += ip->_paramOnLine;
1007             gp_XYZ xyz = lineLoc + ip->_paramOnLine * lineDir;
1008             ip->_node = helper.AddNode( xyz.X(), xyz.Y(), xyz.Z() );
1009             ip->_indexOnLine = nodeCoord-coord0-1;
1010           }
1011           // create a mesh node at ip concident with a grid node
1012           else
1013           {
1014             int nodeIndex = nIndex0 + nShift * ( nodeCoord-coord0 );
1015             if ( !_nodes[ nodeIndex ] )
1016             {
1017               //li.SetIndexOnLine( nodeCoord-coord0 );
1018               //double xyz[3] = { _coords[0][ li.I() ], _coords[1][ li.J() ], _coords[2][ li.K() ]};
1019               gp_XYZ xyz = lineLoc + nodeParam * lineDir;
1020               _nodes   [ nodeIndex ] = helper.AddNode( xyz.X(), xyz.Y(), xyz.Z() );
1021               _gridIntP[ nodeIndex ] = & * ip;
1022             }
1023             if ( _gridIntP[ nodeIndex ] )
1024               _gridIntP[ nodeIndex ]->Add( ip->_faceIDs );
1025             else
1026               _gridIntP[ nodeIndex ] = & * ip;
1027             // ip->_node        = _nodes[ nodeIndex ]; -- to differ from ip on links
1028             ip->_indexOnLine = nodeCoord-coord0;
1029             if ( ++nodeCoord < coordEnd )
1030               nodeParam = *nodeCoord - *coord0;
1031           }
1032         }
1033         // set OUT state to nodes after the last ip
1034         for ( ; nodeCoord < coordEnd; ++nodeCoord )
1035           isNodeOut[ nIndex0 + nShift * ( nodeCoord-coord0 ) ] = true;
1036       }
1037     }
1038
1039     // Create mesh nodes at !OUT nodes of the grid
1040
1041     for ( size_t z = 0; z < _coords[2].size(); ++z )
1042       for ( size_t y = 0; y < _coords[1].size(); ++y )
1043         for ( size_t x = 0; x < _coords[0].size(); ++x )
1044         {
1045           size_t nodeIndex = NodeIndex( x, y, z );
1046           if ( !isNodeOut[ nodeIndex ] && !_nodes[ nodeIndex] )
1047           {
1048             //_nodes[ nodeIndex ] = helper.AddNode( _coords[0][x], _coords[1][y], _coords[2][z] );
1049             gp_XYZ xyz = ( _coords[0][x] * _axes[0] +
1050                            _coords[1][y] * _axes[1] +
1051                            _coords[2][z] * _axes[2] );
1052             _nodes[ nodeIndex ] = helper.AddNode( xyz.X(), xyz.Y(), xyz.Z() );
1053           }
1054         }
1055
1056 #ifdef _MY_DEBUG_
1057     // check validity of transitions
1058     const char* trName[] = { "TANGENT", "IN", "OUT", "APEX" };
1059     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
1060     {
1061       LineIndexer li = GetLineIndexer( iDir );
1062       for ( ; li.More(); ++li )
1063       {
1064         multiset< F_IntersectPoint >& intPnts = _lines[ iDir ][ li.LineIndex() ]._intPoints;
1065         if ( intPnts.empty() ) continue;
1066         if ( intPnts.size() == 1 )
1067         {
1068           if ( intPnts.begin()->_transition != Trans_TANGENT &&
1069                intPnts.begin()->_transition != Trans_APEX )
1070           throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
1071                                     SMESH_Comment("Wrong SOLE transition of GridLine (")
1072                                     << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
1073                                     << ") along " << li._nameConst
1074                                     << ": " << trName[ intPnts.begin()->_transition] );
1075         }
1076         else
1077         {
1078           if ( intPnts.begin()->_transition == Trans_OUT )
1079             throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
1080                                       SMESH_Comment("Wrong START transition of GridLine (")
1081                                       << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
1082                                       << ") along " << li._nameConst
1083                                       << ": " << trName[ intPnts.begin()->_transition ]);
1084           if ( intPnts.rbegin()->_transition == Trans_IN )
1085             throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
1086                                       SMESH_Comment("Wrong END transition of GridLine (")
1087                                       << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
1088                                       << ") along " << li._nameConst
1089                                     << ": " << trName[ intPnts.rbegin()->_transition ]);
1090         }
1091       }
1092     }
1093 #endif
1094   }
1095
1096   //=============================================================================
1097   /*
1098    * Intersects TopoDS_Face with all GridLine's
1099    */
1100   void FaceGridIntersector::Intersect()
1101   {
1102     FaceLineIntersector intersector;
1103     intersector._surfaceInt = GetCurveFaceIntersector();
1104     intersector._tol        = _grid->_tol;
1105     intersector._transOut   = _face.Orientation() == TopAbs_REVERSED ? Trans_IN : Trans_OUT;
1106     intersector._transIn    = _face.Orientation() == TopAbs_REVERSED ? Trans_OUT : Trans_IN;
1107
1108     typedef void (FaceLineIntersector::* PIntFun )(const GridLine& gridLine);
1109     PIntFun interFunction;
1110
1111     bool isDirect = true;
1112     BRepAdaptor_Surface surf( _face );
1113     switch ( surf.GetType() ) {
1114     case GeomAbs_Plane:
1115       intersector._plane = surf.Plane();
1116       interFunction = &FaceLineIntersector::IntersectWithPlane;
1117       isDirect = intersector._plane.Direct();
1118       break;
1119     case GeomAbs_Cylinder:
1120       intersector._cylinder = surf.Cylinder();
1121       interFunction = &FaceLineIntersector::IntersectWithCylinder;
1122       isDirect = intersector._cylinder.Direct();
1123       break;
1124     case GeomAbs_Cone:
1125       intersector._cone = surf.Cone();
1126       interFunction = &FaceLineIntersector::IntersectWithCone;
1127       //isDirect = intersector._cone.Direct();
1128       break;
1129     case GeomAbs_Sphere:
1130       intersector._sphere = surf.Sphere();
1131       interFunction = &FaceLineIntersector::IntersectWithSphere;
1132       isDirect = intersector._sphere.Direct();
1133       break;
1134     case GeomAbs_Torus:
1135       intersector._torus = surf.Torus();
1136       interFunction = &FaceLineIntersector::IntersectWithTorus;
1137       //isDirect = intersector._torus.Direct();
1138       break;
1139     default:
1140       interFunction = &FaceLineIntersector::IntersectWithSurface;
1141     }
1142     if ( !isDirect )
1143       std::swap( intersector._transOut, intersector._transIn );
1144
1145     _intersections.clear();
1146     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
1147     {
1148       if ( surf.GetType() == GeomAbs_Plane )
1149       {
1150         // check if all lines in this direction are parallel to a plane
1151         if ( intersector._plane.Axis().IsNormal( _grid->_lines[iDir][0]._line.Position(),
1152                                                  Precision::Angular()))
1153           continue;
1154         // find out a transition, that is the same for all lines of a direction
1155         gp_Dir plnNorm = intersector._plane.Axis().Direction();
1156         gp_Dir lineDir = _grid->_lines[iDir][0]._line.Direction();
1157         intersector._transition =
1158           ( plnNorm * lineDir < 0 ) ? intersector._transIn : intersector._transOut;
1159       }
1160       if ( surf.GetType() == GeomAbs_Cylinder )
1161       {
1162         // check if all lines in this direction are parallel to a cylinder
1163         if ( intersector._cylinder.Axis().IsParallel( _grid->_lines[iDir][0]._line.Position(),
1164                                                       Precision::Angular()))
1165           continue;
1166       }
1167
1168       // intersect the grid lines with the face
1169       for ( size_t iL = 0; iL < _grid->_lines[iDir].size(); ++iL )
1170       {
1171         GridLine& gridLine = _grid->_lines[iDir][iL];
1172         if ( _bndBox.IsOut( gridLine._line )) continue;
1173
1174         intersector._intPoints.clear();
1175         (intersector.*interFunction)( gridLine ); // <- intersection with gridLine
1176         for ( size_t i = 0; i < intersector._intPoints.size(); ++i )
1177           _intersections.push_back( make_pair( &gridLine, intersector._intPoints[i] ));
1178       }
1179     }
1180   }
1181   //================================================================================
1182   /*
1183    * Return true if (_u,_v) is on the face
1184    */
1185   bool FaceLineIntersector::UVIsOnFace() const
1186   {
1187     TopAbs_State state = _surfaceInt->ClassifyUVPoint(gp_Pnt2d( _u,_v ));
1188     return ( state == TopAbs_IN || state == TopAbs_ON );
1189   }
1190   //================================================================================
1191   /*
1192    * Store an intersection if it is IN or ON the face
1193    */
1194   void FaceLineIntersector::addIntPoint(const bool toClassify)
1195   {
1196     if ( !toClassify || UVIsOnFace() )
1197     {
1198       F_IntersectPoint p;
1199       p._paramOnLine = _w;
1200       p._transition  = _transition;
1201       _intPoints.push_back( p );
1202     }
1203   }
1204   //================================================================================
1205   /*
1206    * Intersect a line with a plane
1207    */
1208   void FaceLineIntersector::IntersectWithPlane   (const GridLine& gridLine)
1209   {
1210     IntAna_IntConicQuad linPlane( gridLine._line, _plane, Precision::Angular());
1211     _w = linPlane.ParamOnConic(1);
1212     if ( isParamOnLineOK( gridLine._length ))
1213     {
1214       ElSLib::Parameters(_plane, linPlane.Point(1) ,_u,_v);
1215       addIntPoint();
1216     }
1217   }
1218   //================================================================================
1219   /*
1220    * Intersect a line with a cylinder
1221    */
1222   void FaceLineIntersector::IntersectWithCylinder(const GridLine& gridLine)
1223   {
1224     IntAna_IntConicQuad linCylinder( gridLine._line, _cylinder );
1225     if ( linCylinder.IsDone() && linCylinder.NbPoints() > 0 )
1226     {
1227       _w = linCylinder.ParamOnConic(1);
1228       if ( linCylinder.NbPoints() == 1 )
1229         _transition = Trans_TANGENT;
1230       else
1231         _transition = _w < linCylinder.ParamOnConic(2) ? _transIn : _transOut;
1232       if ( isParamOnLineOK( gridLine._length ))
1233       {
1234         ElSLib::Parameters(_cylinder, linCylinder.Point(1) ,_u,_v);
1235         addIntPoint();
1236       }
1237       if ( linCylinder.NbPoints() > 1 )
1238       {
1239         _w = linCylinder.ParamOnConic(2);
1240         if ( isParamOnLineOK( gridLine._length ))
1241         {
1242           ElSLib::Parameters(_cylinder, linCylinder.Point(2) ,_u,_v);
1243           _transition = ( _transition == Trans_OUT ) ? Trans_IN : Trans_OUT;
1244           addIntPoint();
1245         }
1246       }
1247     }
1248   }
1249   //================================================================================
1250   /*
1251    * Intersect a line with a cone
1252    */
1253   void FaceLineIntersector::IntersectWithCone (const GridLine& gridLine)
1254   {
1255     IntAna_IntConicQuad linCone(gridLine._line,_cone);
1256     if ( !linCone.IsDone() ) return;
1257     gp_Pnt P;
1258     gp_Vec du, dv, norm;
1259     for ( int i = 1; i <= linCone.NbPoints(); ++i )
1260     {
1261       _w = linCone.ParamOnConic( i );
1262       if ( !isParamOnLineOK( gridLine._length )) continue;
1263       ElSLib::Parameters(_cone, linCone.Point(i) ,_u,_v);
1264       if ( UVIsOnFace() )
1265       {
1266         ElSLib::D1( _u, _v, _cone, P, du, dv );
1267         norm = du ^ dv;
1268         double normSize2 = norm.SquareMagnitude();
1269         if ( normSize2 > Precision::Angular() * Precision::Angular() )
1270         {
1271           double cos = norm.XYZ() * gridLine._line.Direction().XYZ();
1272           cos /= sqrt( normSize2 );
1273           if ( cos < -Precision::Angular() )
1274             _transition = _transIn;
1275           else if ( cos > Precision::Angular() )
1276             _transition = _transOut;
1277           else
1278             _transition = Trans_TANGENT;
1279         }
1280         else
1281         {
1282           _transition = Trans_APEX;
1283         }
1284         addIntPoint( /*toClassify=*/false);
1285       }
1286     }
1287   }
1288   //================================================================================
1289   /*
1290    * Intersect a line with a sphere
1291    */
1292   void FaceLineIntersector::IntersectWithSphere  (const GridLine& gridLine)
1293   {
1294     IntAna_IntConicQuad linSphere(gridLine._line,_sphere);
1295     if ( linSphere.IsDone() && linSphere.NbPoints() > 0 )
1296     {
1297       _w = linSphere.ParamOnConic(1);
1298       if ( linSphere.NbPoints() == 1 )
1299         _transition = Trans_TANGENT;
1300       else
1301         _transition = _w < linSphere.ParamOnConic(2) ? _transIn : _transOut;
1302       if ( isParamOnLineOK( gridLine._length ))
1303       {
1304         ElSLib::Parameters(_sphere, linSphere.Point(1) ,_u,_v);
1305         addIntPoint();
1306       }
1307       if ( linSphere.NbPoints() > 1 )
1308       {
1309         _w = linSphere.ParamOnConic(2);
1310         if ( isParamOnLineOK( gridLine._length ))
1311         {
1312           ElSLib::Parameters(_sphere, linSphere.Point(2) ,_u,_v);
1313           _transition = ( _transition == Trans_OUT ) ? Trans_IN : Trans_OUT;
1314           addIntPoint();
1315         }
1316       }
1317     }
1318   }
1319   //================================================================================
1320   /*
1321    * Intersect a line with a torus
1322    */
1323   void FaceLineIntersector::IntersectWithTorus   (const GridLine& gridLine)
1324   {
1325     IntAna_IntLinTorus linTorus(gridLine._line,_torus);
1326     if ( !linTorus.IsDone()) return;
1327     gp_Pnt P;
1328     gp_Vec du, dv, norm;
1329     for ( int i = 1; i <= linTorus.NbPoints(); ++i )
1330     {
1331       _w = linTorus.ParamOnLine( i );
1332       if ( !isParamOnLineOK( gridLine._length )) continue;
1333       linTorus.ParamOnTorus( i, _u,_v );
1334       if ( UVIsOnFace() )
1335       {
1336         ElSLib::D1( _u, _v, _torus, P, du, dv );
1337         norm = du ^ dv;
1338         double normSize = norm.Magnitude();
1339         double cos = norm.XYZ() * gridLine._line.Direction().XYZ();
1340         cos /= normSize;
1341         if ( cos < -Precision::Angular() )
1342           _transition = _transIn;
1343         else if ( cos > Precision::Angular() )
1344           _transition = _transOut;
1345         else
1346           _transition = Trans_TANGENT;
1347         addIntPoint( /*toClassify=*/false);
1348       }
1349     }
1350   }
1351   //================================================================================
1352   /*
1353    * Intersect a line with a non-analytical surface
1354    */
1355   void FaceLineIntersector::IntersectWithSurface (const GridLine& gridLine)
1356   {
1357     _surfaceInt->Perform( gridLine._line, 0.0, gridLine._length );
1358     if ( !_surfaceInt->IsDone() ) return;
1359     for ( int i = 1; i <= _surfaceInt->NbPnt(); ++i )
1360     {
1361       _transition = Transition( _surfaceInt->Transition( i ) );
1362       _w = _surfaceInt->WParameter( i );
1363       addIntPoint(/*toClassify=*/false);
1364     }
1365   }
1366   //================================================================================
1367   /*
1368    * check if its face can be safely intersected in a thread
1369    */
1370   bool FaceGridIntersector::IsThreadSafe(set< const Standard_Transient* >& noSafeTShapes) const
1371   {
1372     bool isSafe = true;
1373
1374     // check surface
1375     TopLoc_Location loc;
1376     Handle(Geom_Surface) surf = BRep_Tool::Surface( _face, loc );
1377     Handle(Geom_RectangularTrimmedSurface) ts =
1378       Handle(Geom_RectangularTrimmedSurface)::DownCast( surf );
1379     while( !ts.IsNull() ) {
1380       surf = ts->BasisSurface();
1381       ts = Handle(Geom_RectangularTrimmedSurface)::DownCast(surf);
1382     }
1383     if ( surf->IsKind( STANDARD_TYPE(Geom_BSplineSurface )) ||
1384          surf->IsKind( STANDARD_TYPE(Geom_BezierSurface )))
1385       if ( !noSafeTShapes.insert((const Standard_Transient*) _face.TShape() ).second )
1386         isSafe = false;
1387
1388     double f, l;
1389     TopExp_Explorer exp( _face, TopAbs_EDGE );
1390     for ( ; exp.More(); exp.Next() )
1391     {
1392       bool edgeIsSafe = true;
1393       const TopoDS_Edge& e = TopoDS::Edge( exp.Current() );
1394       // check 3d curve
1395       {
1396         Handle(Geom_Curve) c = BRep_Tool::Curve( e, loc, f, l);
1397         if ( !c.IsNull() )
1398         {
1399           Handle(Geom_TrimmedCurve) tc = Handle(Geom_TrimmedCurve)::DownCast(c);
1400           while( !tc.IsNull() ) {
1401             c = tc->BasisCurve();
1402             tc = Handle(Geom_TrimmedCurve)::DownCast(c);
1403           }
1404           if ( c->IsKind( STANDARD_TYPE(Geom_BSplineCurve )) ||
1405                c->IsKind( STANDARD_TYPE(Geom_BezierCurve )))
1406             edgeIsSafe = false;
1407         }
1408       }
1409       // check 2d curve
1410       if ( edgeIsSafe )
1411       {
1412         Handle(Geom2d_Curve) c2 = BRep_Tool::CurveOnSurface( e, surf, loc, f, l);
1413         if ( !c2.IsNull() )
1414         {
1415           Handle(Geom2d_TrimmedCurve) tc = Handle(Geom2d_TrimmedCurve)::DownCast(c2);
1416           while( !tc.IsNull() ) {
1417             c2 = tc->BasisCurve();
1418             tc = Handle(Geom2d_TrimmedCurve)::DownCast(c2);
1419           }
1420           if ( c2->IsKind( STANDARD_TYPE(Geom2d_BSplineCurve )) ||
1421                c2->IsKind( STANDARD_TYPE(Geom2d_BezierCurve )))
1422             edgeIsSafe = false;
1423         }
1424       }
1425       if ( !edgeIsSafe && !noSafeTShapes.insert((const Standard_Transient*) e.TShape() ).second )
1426         isSafe = false;
1427     }
1428     return isSafe;
1429   }
1430   //================================================================================
1431   /*!
1432    * \brief Creates topology of the hexahedron
1433    */
1434   Hexahedron::Hexahedron(const double sizeThreshold, Grid* grid)
1435     : _grid( grid ), _sizeThreshold( sizeThreshold ), _nbFaceIntNodes(0)
1436   {
1437     _polygons.reserve(100); // to avoid reallocation;
1438
1439     //set nodes shift within grid->_nodes from the node 000 
1440     size_t dx = _grid->NodeIndexDX();
1441     size_t dy = _grid->NodeIndexDY();
1442     size_t dz = _grid->NodeIndexDZ();
1443     size_t i000 = 0;
1444     size_t i100 = i000 + dx;
1445     size_t i010 = i000 + dy;
1446     size_t i110 = i010 + dx;
1447     size_t i001 = i000 + dz;
1448     size_t i101 = i100 + dz;
1449     size_t i011 = i010 + dz;
1450     size_t i111 = i110 + dz;
1451     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V000 )] = i000;
1452     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V100 )] = i100;
1453     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V010 )] = i010;
1454     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V110 )] = i110;
1455     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V001 )] = i001;
1456     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V101 )] = i101;
1457     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V011 )] = i011;
1458     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V111 )] = i111;
1459
1460     vector< int > idVec;
1461     // set nodes to links
1462     for ( int linkID = SMESH_Block::ID_Ex00; linkID <= SMESH_Block::ID_E11z; ++linkID )
1463     {
1464       SMESH_Block::GetEdgeVertexIDs( linkID, idVec );
1465       _Link& link = _hexLinks[ SMESH_Block::ShapeIndex( linkID )];
1466       link._nodes[0] = &_hexNodes[ SMESH_Block::ShapeIndex( idVec[0] )];
1467       link._nodes[1] = &_hexNodes[ SMESH_Block::ShapeIndex( idVec[1] )];
1468     }
1469
1470     // set links to faces
1471     int interlace[4] = { 0, 3, 1, 2 }; // to walk by links around a face: { u0, 1v, u1, 0v }
1472     for ( int faceID = SMESH_Block::ID_Fxy0; faceID <= SMESH_Block::ID_F1yz; ++faceID )
1473     {
1474       SMESH_Block::GetFaceEdgesIDs( faceID, idVec );
1475       _Face& quad = _hexQuads[ SMESH_Block::ShapeIndex( faceID )];
1476       bool revFace = ( faceID == SMESH_Block::ID_Fxy0 ||
1477                        faceID == SMESH_Block::ID_Fx1z ||
1478                        faceID == SMESH_Block::ID_F0yz );
1479       quad._links.resize(4);
1480       vector<_OrientedLink>::iterator         frwLinkIt = quad._links.begin();
1481       vector<_OrientedLink>::reverse_iterator revLinkIt = quad._links.rbegin();
1482       for ( int i = 0; i < 4; ++i )
1483       {
1484         bool revLink = revFace;
1485         if ( i > 1 ) // reverse links u1 and v0
1486           revLink = !revLink;
1487         _OrientedLink& link = revFace ? *revLinkIt++ : *frwLinkIt++;
1488         link = _OrientedLink( & _hexLinks[ SMESH_Block::ShapeIndex( idVec[interlace[i]] )],
1489                               revLink );
1490       }
1491     }
1492   }
1493   //================================================================================
1494   /*!
1495    * \brief Copy constructor
1496    */
1497   Hexahedron::Hexahedron( const Hexahedron& other )
1498     :_grid( other._grid ), _sizeThreshold( other._sizeThreshold ), _nbFaceIntNodes(0)
1499   {
1500     _polygons.reserve(100); // to avoid reallocation;
1501
1502     for ( int i = 0; i < 8; ++i )
1503       _nodeShift[i] = other._nodeShift[i];
1504
1505     for ( int i = 0; i < 12; ++i )
1506     {
1507       const _Link& srcLink = other._hexLinks[ i ];
1508       _Link&       tgtLink = this->_hexLinks[ i ];
1509       tgtLink._nodes[0] = _hexNodes + ( srcLink._nodes[0] - other._hexNodes );
1510       tgtLink._nodes[1] = _hexNodes + ( srcLink._nodes[1] - other._hexNodes );
1511     }
1512
1513     for ( int i = 0; i < 6; ++i )
1514     {
1515       const _Face& srcQuad = other._hexQuads[ i ];
1516       _Face&       tgtQuad = this->_hexQuads[ i ];
1517       tgtQuad._links.resize(4);
1518       for ( int j = 0; j < 4; ++j )
1519       {
1520         const _OrientedLink& srcLink = srcQuad._links[ j ];
1521         _OrientedLink&       tgtLink = tgtQuad._links[ j ];
1522         tgtLink._reverse = srcLink._reverse;
1523         tgtLink._link    = _hexLinks + ( srcLink._link - other._hexLinks );
1524       }
1525     }
1526   }
1527   
1528   //================================================================================
1529   /*!
1530    * \brief Initializes its data by given grid cell
1531    */
1532   void Hexahedron::init( size_t i, size_t j, size_t k )
1533   {
1534     _i = i; _j = j; _k = k;
1535     // set nodes of grid to nodes of the hexahedron and
1536     // count nodes at hexahedron corners located IN and ON geometry
1537     _nbCornerNodes = _nbBndNodes = 0;
1538     _origNodeInd   = _grid->NodeIndex( i,j,k );
1539     for ( int iN = 0; iN < 8; ++iN )
1540     {
1541       _hexNodes[iN]._node     = _grid->_nodes   [ _origNodeInd + _nodeShift[iN] ];
1542       _hexNodes[iN]._intPoint = _grid->_gridIntP[ _origNodeInd + _nodeShift[iN] ];
1543       _nbCornerNodes += bool( _hexNodes[iN]._node );
1544       _nbBndNodes    += bool( _hexNodes[iN]._intPoint );
1545     }
1546
1547     _sideLength[0] = _grid->_coords[0][i+1] - _grid->_coords[0][i];
1548     _sideLength[1] = _grid->_coords[1][j+1] - _grid->_coords[1][j];
1549     _sideLength[2] = _grid->_coords[2][k+1] - _grid->_coords[2][k];
1550
1551     _intNodes.clear();
1552     _vIntNodes.clear();
1553
1554     if ( _nbFaceIntNodes + _eIntPoints.size() > 0 &&
1555          _nbFaceIntNodes + _nbCornerNodes + _eIntPoints.size() > 3)
1556     {
1557       _intNodes.reserve( 3 * _nbBndNodes + _nbFaceIntNodes + _eIntPoints.size() );
1558
1559       _Link split;
1560       // create sub-links (_splits) by splitting links with _fIntPoints
1561       for ( int iLink = 0; iLink < 12; ++iLink )
1562       {
1563         _Link& link = _hexLinks[ iLink ];
1564         link._fIntNodes.resize( link._fIntPoints.size() );
1565         for ( size_t i = 0; i < link._fIntPoints.size(); ++i )
1566         {
1567           _intNodes.push_back( _Node( 0, link._fIntPoints[i] ));
1568           link._fIntNodes[ i ] = & _intNodes.back();
1569         }
1570
1571         link._splits.clear();
1572         split._nodes[ 0 ] = link._nodes[0];
1573         bool isOut = ( ! link._nodes[0]->Node() ); // is1stNodeOut( iLink );
1574         bool checkTransition;
1575         for ( size_t i = 0; i < link._fIntNodes.size(); ++i )
1576         {
1577           if ( link._fIntNodes[i]->Node() ) // intersection non-coinsident with a grid node
1578           {
1579             if ( split._nodes[ 0 ]->Node() && !isOut )
1580             {
1581               split._nodes[ 1 ] = link._fIntNodes[i];
1582               link._splits.push_back( split );
1583             }
1584             split._nodes[ 0 ] = link._fIntNodes[i];
1585             checkTransition = true;
1586           }
1587           else // FACE intersection coinsident with a grid node
1588           {
1589             checkTransition = ( link._nodes[0]->Node() );
1590           }
1591           if ( checkTransition )
1592           {
1593             switch ( link._fIntPoints[i]->_transition ) {
1594             case Trans_OUT: isOut = true; break;
1595             case Trans_IN : isOut = false; break;
1596             default:
1597               if ( !link._fIntNodes[i]->Node() && i == 0 )
1598                 isOut = is1stNodeOut( link );
1599               else
1600                 ; // isOut remains the same
1601             }
1602           }
1603         }
1604         if ( link._nodes[ 1 ]->Node() && split._nodes[ 0 ]->Node() && !isOut )
1605         {
1606           split._nodes[ 1 ] = link._nodes[1];
1607           link._splits.push_back( split );
1608         }
1609       }
1610
1611       // Create _Node's at intersections with EDGEs.
1612
1613       const double tol2 = _grid->_tol * _grid->_tol;
1614       int facets[3], nbFacets, subEntity;
1615
1616       for ( size_t iP = 0; iP < _eIntPoints.size(); ++iP )
1617       {
1618         nbFacets = getEntity( _eIntPoints[iP], facets, subEntity );
1619         _Node* equalNode = 0;
1620         switch( nbFacets ) {
1621         case 1: // in a _Face
1622         {
1623           _Face& quad = _hexQuads[ facets[0] - SMESH_Block::ID_FirstF ];
1624           equalNode = FindEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
1625           if ( equalNode ) {
1626             equalNode->Add( _eIntPoints[ iP ] );
1627           }
1628           else {
1629             _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
1630             quad._eIntNodes.push_back( & _intNodes.back() );
1631           }
1632           break;
1633         }
1634         case 2: // on a _Link
1635         {
1636           _Link& link = _hexLinks[ subEntity - SMESH_Block::ID_FirstE ];
1637           if ( link._splits.size() > 0 )
1638           {
1639             equalNode = FindEqualNode( link._fIntNodes, _eIntPoints[ iP ], tol2 );
1640             if ( equalNode )
1641               equalNode->Add( _eIntPoints[ iP ] );
1642           }
1643           else
1644           {
1645             _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
1646             for ( int iF = 0; iF < 2; ++iF )
1647             {
1648               _Face& quad = _hexQuads[ facets[iF] - SMESH_Block::ID_FirstF ];
1649               equalNode = FindEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
1650               if ( equalNode ) {
1651                 equalNode->Add( _eIntPoints[ iP ] );
1652               }
1653               else {
1654                 quad._eIntNodes.push_back( & _intNodes.back() );
1655               }
1656             }
1657           }
1658           break;
1659         }
1660         case 3: // at a corner
1661         {
1662           _Node& node = _hexNodes[ subEntity - SMESH_Block::ID_FirstV ];
1663           if ( node.Node() > 0 )
1664           {
1665             if ( node._intPoint )
1666               node._intPoint->Add( _eIntPoints[ iP ]->_faceIDs, _eIntPoints[ iP ]->_node );
1667           }
1668           else
1669           {
1670             _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
1671             for ( int iF = 0; iF < 3; ++iF )
1672             {
1673               _Face& quad = _hexQuads[ facets[iF] - SMESH_Block::ID_FirstF ];
1674               equalNode = FindEqualNode( quad._eIntNodes, _eIntPoints[ iP ], tol2 );
1675               if ( equalNode ) {
1676                 equalNode->Add( _eIntPoints[ iP ] );
1677               }
1678               else {
1679                 quad._eIntNodes.push_back( & _intNodes.back() );
1680               }
1681             }
1682           }
1683           break;
1684         }
1685         } // switch( nbFacets )
1686
1687         if ( nbFacets == 0 ||
1688              _grid->_shapes( _eIntPoints[ iP ]->_shapeID ).ShapeType() == TopAbs_VERTEX )
1689         {
1690           equalNode = FindEqualNode( _vIntNodes, _eIntPoints[ iP ], tol2 );
1691           if ( equalNode ) {
1692             equalNode->Add( _eIntPoints[ iP ] );
1693           }
1694           else {
1695             if ( _intNodes.empty() || _intNodes.back().EdgeIntPnt() != _eIntPoints[ iP ])
1696               _intNodes.push_back( _Node( 0, _eIntPoints[ iP ]));
1697             _vIntNodes.push_back( & _intNodes.back() );
1698           }
1699         }
1700       } // loop on _eIntPoints
1701     }
1702     else if ( 3 < _nbCornerNodes && _nbCornerNodes < 8 ) // _nbFaceIntNodes == 0
1703     {
1704       _Link split;
1705       // create sub-links (_splits) of whole links
1706       for ( int iLink = 0; iLink < 12; ++iLink )
1707       {
1708         _Link& link = _hexLinks[ iLink ];
1709         link._splits.clear();
1710         if ( link._nodes[ 0 ]->Node() && link._nodes[ 1 ]->Node() )
1711         {
1712           split._nodes[ 0 ] = link._nodes[0];
1713           split._nodes[ 1 ] = link._nodes[1];
1714           link._splits.push_back( split );
1715         }
1716       }
1717     }
1718
1719   }
1720   //================================================================================
1721   /*!
1722    * \brief Initializes its data by given grid cell (countered from zero)
1723    */
1724   void Hexahedron::init( size_t iCell )
1725   {
1726     size_t iNbCell = _grid->_coords[0].size() - 1;
1727     size_t jNbCell = _grid->_coords[1].size() - 1;
1728     _i = iCell % iNbCell;
1729     _j = ( iCell % ( iNbCell * jNbCell )) / iNbCell;
1730     _k = iCell / iNbCell / jNbCell;
1731     init( _i, _j, _k );
1732   }
1733
1734   //================================================================================
1735   /*!
1736    * \brief Compute mesh volumes resulted from intersection of the Hexahedron
1737    */
1738   void Hexahedron::ComputeElements()
1739   {
1740     Init();
1741
1742     int nbIntersections = _nbFaceIntNodes + _eIntPoints.size();
1743     if ( _nbCornerNodes + nbIntersections < 4 )
1744       return;
1745
1746     if ( _nbBndNodes == _nbCornerNodes && nbIntersections == 0 && isInHole() )
1747       return;
1748
1749     _polygons.clear();
1750     _polygons.reserve( 20 );
1751
1752     // Create polygons from quadrangles
1753     // --------------------------------
1754
1755     _Link                   polyLink;
1756     vector< _OrientedLink > splits;
1757     vector<_Node*>          chainNodes, usedEdgeNodes;
1758     _Face*                  coplanarPolyg;
1759
1760     bool hasEdgeIntersections = !_eIntPoints.empty();
1761
1762     for ( int iF = 0; iF < 6; ++iF ) // loop on 6 sides of a hexahedron
1763     {
1764       _Face& quad = _hexQuads[ iF ] ;
1765
1766       _polygons.resize( _polygons.size() + 1 );
1767       _Face* polygon = &_polygons.back();
1768       polygon->_polyLinks.reserve( 20 );
1769
1770       splits.clear();
1771       for ( int iE = 0; iE < 4; ++iE ) // loop on 4 sides of a quadrangle
1772         for ( int iS = 0; iS < quad._links[ iE ].NbResultLinks(); ++iS )
1773           splits.push_back( quad._links[ iE ].ResultLink( iS ));
1774
1775       // add splits of links to a polygon and add _polyLinks to make
1776       // polygon's boundary closed
1777
1778       int nbSplits = splits.size();
1779       if ( nbSplits < 2 && quad._eIntNodes.empty() )
1780         nbSplits = 0;
1781
1782 #ifdef _DEBUG_
1783       for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
1784         if ( quad._eIntNodes[ iP ]->IsUsedInFace( polygon ))
1785           quad._eIntNodes[ iP ]->_usedInFace = 0;
1786 #endif
1787       int nbUsedEdgeNodes = 0;
1788
1789       while ( nbSplits > 0 )
1790       {
1791         size_t iS = 0;
1792         while ( !splits[ iS ] )
1793           ++iS;
1794
1795         if ( !polygon->_links.empty() )
1796         {
1797           _polygons.resize( _polygons.size() + 1 );
1798           polygon = &_polygons.back();
1799           polygon->_polyLinks.reserve( 20 );
1800         }
1801         polygon->_links.push_back( splits[ iS ] );
1802         splits[ iS++ ]._link = 0;
1803         --nbSplits;
1804
1805         _Node* nFirst = polygon->_links.back().FirstNode();
1806         _Node *n1,*n2 = polygon->_links.back().LastNode();
1807         for ( ; nFirst != n2 && iS < splits.size(); ++iS )
1808         {
1809           _OrientedLink& split = splits[ iS ];
1810           if ( !split ) continue;
1811
1812           n1 = split.FirstNode();
1813           if ( n1 != n2 )
1814           {
1815             // try to connect to intersections with EDGEs
1816             if ( quad._eIntNodes.size() > nbUsedEdgeNodes  &&
1817                  findChain( n2, n1, quad, chainNodes ))
1818             {
1819               for ( size_t i = 1; i < chainNodes.size(); ++i )
1820               {
1821                 polyLink._nodes[0] = chainNodes[i-1];
1822                 polyLink._nodes[1] = chainNodes[i];
1823                 polygon->_polyLinks.push_back( polyLink );
1824                 polygon->_links.push_back( _OrientedLink( &polygon->_polyLinks.back() ));
1825                 nbUsedEdgeNodes += ( polyLink._nodes[1]->IsUsedInFace( polygon ));
1826               }
1827               if ( chainNodes.back() != n1 )
1828               {
1829                 n2 = chainNodes.back();
1830                 --iS;
1831                 continue;
1832               }
1833             }
1834             // try to connect to a split ending on the same FACE
1835             else
1836             {
1837               _OrientedLink foundSplit;
1838               for ( int i = iS; i < splits.size() && !foundSplit; ++i )
1839                 if (( foundSplit = splits[ i ]) &&
1840                     ( n2->IsLinked( foundSplit.FirstNode()->_intPoint )))
1841                 {
1842                   polyLink._nodes[0] = n2;
1843                   polyLink._nodes[1] = foundSplit.FirstNode();
1844                   polygon->_polyLinks.push_back( polyLink );
1845                   polygon->_links.push_back( _OrientedLink( &polygon->_polyLinks.back() ));
1846                   iS = i - 1;
1847                 }
1848                 else
1849                 {
1850                   foundSplit._link = 0;
1851                 }
1852               if ( foundSplit )
1853               {
1854                 n2 = foundSplit.FirstNode();
1855                 continue;
1856               }
1857               else
1858               {
1859                 if ( n2->IsLinked( nFirst->_intPoint ))
1860                   break;
1861                 polyLink._nodes[0] = n2;
1862                 polyLink._nodes[1] = n1;
1863                 polygon->_polyLinks.push_back( polyLink );
1864                 polygon->_links.push_back( _OrientedLink( &polygon->_polyLinks.back() ));
1865               }
1866             }
1867           }
1868           polygon->_links.push_back( split );
1869           split._link = 0;
1870           --nbSplits;
1871           n2 = polygon->_links.back().LastNode();
1872
1873         } // loop on splits
1874
1875         if ( nFirst != n2 ) // close a polygon
1876         {
1877           if ( !findChain( n2, nFirst, quad, chainNodes ))
1878           {
1879             if ( !closePolygon( polygon, chainNodes ))
1880               chainNodes.push_back( nFirst );
1881           }
1882           for ( size_t i = 1; i < chainNodes.size(); ++i )
1883           {
1884             polyLink._nodes[0] = chainNodes[i-1];
1885             polyLink._nodes[1] = chainNodes[i];
1886             polygon->_polyLinks.push_back( polyLink );
1887             polygon->_links.push_back( _OrientedLink( &polygon->_polyLinks.back() ));
1888             nbUsedEdgeNodes += bool( polyLink._nodes[1]->IsUsedInFace( polygon ));
1889           }
1890         }
1891
1892         if ( polygon->_links.size() < 3 && nbSplits > 0 )
1893         {
1894           polygon->_polyLinks.clear();
1895           polygon->_links.clear();
1896         }
1897       } // while ( nbSplits > 0 )
1898
1899       // if ( quad._eIntNodes.size() > nbUsedEdgeNodes )
1900       // {
1901       //   // make _vIntNodes from not used _eIntNodes
1902       //   const double tol = 0.05 * Min( Min( _sideLength[0], _sideLength[1] ), _sideLength[0] );
1903       //   for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
1904       //   {
1905       //     if ( quad._eIntNodes[ iP ]->IsUsedInFace() ) continue;
1906       //     _Node* equalNode =
1907       //       FindEqualNode( _vIntNodes, quad._eIntNodes[ iP ].EdgeIntPnt(), tol*tol );
1908       //     if ( equalNode )
1909       //       equalNode->Add( quad._eIntNodes[ iP ].EdgeIntPnt() );
1910       //     else
1911       //       _vIntNodes.push_back( quad._eIntNodes[ iP ]);
1912       //   }
1913       // }
1914
1915       if ( polygon->_links.size() < 3 )
1916       {
1917         _polygons.pop_back();
1918         //usedEdgeNodes.resize( usedEdgeNodes.size() - nbUsedEdgeNodes );
1919       }
1920     }  // loop on 6 hexahedron sides
1921
1922     // Create polygons closing holes in a polyhedron
1923     // ----------------------------------------------
1924
1925     // clear _usedInFace
1926     for ( size_t iN = 0; iN < _intNodes.size(); ++iN )
1927       _intNodes[ iN ]._usedInFace = 0;
1928
1929     // add polygons to their links and mark used nodes
1930     for ( size_t iP = 0; iP < _polygons.size(); ++iP )
1931     {
1932       _Face& polygon = _polygons[ iP ];
1933       for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
1934       {
1935         polygon._links[ iL ].AddFace( &polygon );
1936         polygon._links[ iL ].FirstNode()->_usedInFace = &polygon;
1937       }
1938     }
1939     // find free links
1940     vector< _OrientedLink* > freeLinks;
1941     freeLinks.reserve(20);
1942     for ( size_t iP = 0; iP < _polygons.size(); ++iP )
1943     {
1944       _Face& polygon = _polygons[ iP ];
1945       for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
1946         if ( polygon._links[ iL ].NbFaces() < 2 )
1947         {
1948           freeLinks.push_back( & polygon._links[ iL ]);
1949           freeLinks.back()->FirstNode()->IsUsedInFace() == true;
1950         }
1951     }
1952     int nbFreeLinks = freeLinks.size();
1953     if ( nbFreeLinks > 0 && nbFreeLinks < 3 ) return;
1954
1955     // put not used intersection nodes to _vIntNodes
1956     int nbVertexNodes = 0; // nb not used vertex nodes
1957     {
1958       for ( size_t iN = 0; iN < _vIntNodes.size(); ++iN )
1959         nbVertexNodes += ( !_vIntNodes[ iN ]->IsUsedInFace() );
1960
1961       const double tol = 1e-3 * Min( Min( _sideLength[0], _sideLength[1] ), _sideLength[0] );
1962       for ( size_t iN = _nbFaceIntNodes; iN < _intNodes.size(); ++iN )
1963       {
1964         if ( _intNodes[ iN ].IsUsedInFace() ) continue;
1965         if ( dynamic_cast< const F_IntersectPoint* >( _intNodes[ iN ]._intPoint )) continue;
1966         _Node* equalNode =
1967           FindEqualNode( _vIntNodes, _intNodes[ iN ].EdgeIntPnt(), tol*tol );
1968         if ( !equalNode /*|| equalNode->IsUsedInFace()*/ )
1969         {
1970           _vIntNodes.push_back( &_intNodes[ iN ]);
1971           ++nbVertexNodes;
1972         }
1973       }
1974     }
1975
1976     set<TGeomID> usedFaceIDs;
1977     TGeomID curFace = 0;
1978     const size_t nbQuadPolygons = _polygons.size();
1979
1980     // create polygons by making closed chains of free links
1981     size_t iPolygon = _polygons.size();
1982     while ( nbFreeLinks > 0 )
1983     {
1984       if ( iPolygon == _polygons.size() )
1985         _polygons.resize( _polygons.size() + 1 );
1986       _Face& polygon = _polygons[ iPolygon ];
1987       polygon._polyLinks.reserve( 20 );
1988       polygon._links.reserve( 20 );
1989
1990       _OrientedLink* curLink = 0;
1991       _Node*         curNode;
1992       if (( !hasEdgeIntersections ) ||
1993           ( nbFreeLinks < 4 && nbVertexNodes == 0 ))
1994       {
1995         // get a remaining link to start from
1996         for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
1997           if (( curLink = freeLinks[ iL ] ))
1998             freeLinks[ iL ] = 0;
1999         polygon._links.push_back( *curLink );
2000         --nbFreeLinks;
2001         do
2002         {
2003           // find all links connected to curLink
2004           curNode = curLink->FirstNode();
2005           curLink = 0;
2006           for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2007             if ( freeLinks[ iL ] && freeLinks[ iL ]->LastNode() == curNode )
2008             {
2009               curLink = freeLinks[ iL ];
2010               freeLinks[ iL ] = 0;
2011               --nbFreeLinks;
2012               polygon._links.push_back( *curLink );
2013             }
2014         } while ( curLink );
2015       }
2016       else // there are intersections with EDGEs
2017       {
2018         // get a remaining link to start from, one lying on minimal
2019         // nb of FACEs
2020         {
2021           vector< pair< TGeomID, int > > facesOfLink[3];
2022           pair< TGeomID, int > faceOfLink( -1, -1 );
2023           vector< TGeomID > faces;
2024           for ( size_t iL = 0; iL < freeLinks.size(); ++iL )
2025             if ( freeLinks[ iL ] )
2026             {
2027               faces = freeLinks[ iL ]->GetNotUsedFace( usedFaceIDs );
2028               if ( faces.size() == 1 )
2029               {
2030                 faceOfLink = make_pair( faces[0], iL );
2031                 if ( !freeLinks[ iL ]->HasEdgeNodes() )
2032                   break;
2033                 facesOfLink[0].push_back( faceOfLink );
2034               }
2035               else if ( facesOfLink[0].empty() )
2036               {
2037                 faceOfLink = make_pair(( faces.empty() ? -1 : faces[0]), iL );
2038                 facesOfLink[ 1 + faces.empty() ].push_back( faceOfLink );
2039               }
2040             }
2041           for ( int i = 0; faceOfLink.second < 0 && i < 3; ++i )
2042             if ( !facesOfLink[i].empty() )
2043               faceOfLink = facesOfLink[i][0];
2044
2045           if ( faceOfLink.first < 0 ) // all faces used
2046           {
2047             for ( size_t i = 0; i < facesOfLink[2].size() && faceOfLink.first < 1; ++i )
2048             {
2049               curLink = freeLinks[ facesOfLink[2][i].second ];
2050               faceOfLink.first = curLink->FirstNode()->IsLinked( curLink->LastNode()->_intPoint );
2051             }
2052             usedFaceIDs.clear();
2053           }
2054           curFace = faceOfLink.first;
2055           curLink = freeLinks[ faceOfLink.second ];
2056           freeLinks[ faceOfLink.second ] = 0;
2057         }
2058         usedFaceIDs.insert( curFace );
2059         polygon._links.push_back( *curLink );
2060         --nbFreeLinks;
2061
2062         // find all links bounding a FACE of curLink
2063         do
2064         {
2065           // go forward from curLink
2066           curNode = curLink->LastNode();
2067           curLink = 0;
2068           for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2069             if ( freeLinks[ iL ] &&
2070                  freeLinks[ iL ]->FirstNode() == curNode &&
2071                  freeLinks[ iL ]->LastNode()->IsOnFace( curFace ))
2072             {
2073               curLink = freeLinks[ iL ];
2074               freeLinks[ iL ] = 0;
2075               polygon._links.push_back( *curLink );
2076               --nbFreeLinks;
2077             }
2078         } while ( curLink );
2079
2080         std::reverse( polygon._links.begin(), polygon._links.end() );
2081
2082         curLink = & polygon._links.back();
2083         do
2084         {
2085           // go backward from curLink
2086           curNode = curLink->FirstNode();
2087           curLink = 0;
2088           for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
2089             if ( freeLinks[ iL ] &&
2090                  freeLinks[ iL ]->LastNode() == curNode &&
2091                  freeLinks[ iL ]->FirstNode()->IsOnFace( curFace ))
2092             {
2093               curLink = freeLinks[ iL ];
2094               freeLinks[ iL ] = 0;
2095               polygon._links.push_back( *curLink );
2096               --nbFreeLinks;
2097             }
2098         } while ( curLink );
2099
2100         curNode = polygon._links.back().FirstNode();
2101
2102         if ( polygon._links[0].LastNode() != curNode )
2103         {
2104           if ( nbVertexNodes > 0 )
2105           {
2106             // add links with _vIntNodes if not already used
2107             for ( size_t iN = 0; iN < _vIntNodes.size(); ++iN )
2108               if ( !_vIntNodes[ iN ]->IsUsedInFace() &&
2109                    _vIntNodes[ iN ]->IsOnFace( curFace ))
2110               {
2111                 _vIntNodes[ iN ]->_usedInFace = &polygon;
2112                 --nbVertexNodes;
2113                 polyLink._nodes[0] = _vIntNodes[ iN ];
2114                 polyLink._nodes[1] = curNode;
2115                 polygon._polyLinks.push_back( polyLink );
2116                 polygon._links.push_back( _OrientedLink( &polygon._polyLinks.back() ));
2117                 freeLinks.push_back( &polygon._links.back() );
2118                 ++nbFreeLinks;
2119                 curNode = _vIntNodes[ iN ];
2120                 // TODO: to reorder _vIntNodes within polygon, if there are several ones
2121               }
2122           }
2123           // if ( polygon._links.size() > 1 )
2124           {
2125             polyLink._nodes[0] = polygon._links[0].LastNode();
2126             polyLink._nodes[1] = curNode;
2127             polygon._polyLinks.push_back( polyLink );
2128             polygon._links.push_back( _OrientedLink( &polygon._polyLinks.back() ));
2129             freeLinks.push_back( &polygon._links.back() );
2130             ++nbFreeLinks;
2131           }
2132         }
2133       } // if there are intersections with EDGEs
2134
2135       if ( polygon._links.size() < 2 ||
2136            polygon._links[0].LastNode() != polygon._links.back().FirstNode() )
2137         return; // closed polygon not found -> invalid polyhedron
2138
2139       if ( polygon._links.size() == 2 )
2140       {
2141         if ( freeLinks.back() == &polygon._links.back() )
2142         {
2143           freeLinks.pop_back();
2144           --nbFreeLinks;
2145         }
2146         if ( polygon._links.front().NbFaces() > 0 )
2147           polygon._links.back().AddFace( polygon._links.front()._link->_faces[0] );
2148         if ( polygon._links.back().NbFaces() > 0 )
2149           polygon._links.front().AddFace( polygon._links.back()._link->_faces[0] );
2150
2151         _polygons.pop_back();
2152       }
2153       else // polygon._links.size() >= 2
2154       {
2155         // add polygon to its links
2156         for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
2157         {
2158           polygon._links[ iL ].AddFace( &polygon );
2159           polygon._links[ iL ].Reverse();
2160         }
2161         if ( hasEdgeIntersections && iPolygon == _polygons.size() - 1 )
2162         {
2163           // check that a polygon does not lie in the plane of another polygon
2164           coplanarPolyg = 0;
2165           for ( size_t iL = 0; iL < polygon._links.size() && !coplanarPolyg; ++iL )
2166           {
2167             if ( polygon._links[ iL ].NbFaces() < 2 )
2168               continue; // it's a just added free link
2169             // look for a polygon made on a hexa side and sharing
2170             // two or more haxa links
2171             size_t iL2;
2172             coplanarPolyg = polygon._links[ iL ]._link->_faces[0];
2173             for ( iL2 = iL + 1; iL2 < polygon._links.size(); ++iL2 )
2174               if ( polygon._links[ iL2 ]._link->_faces[0] == coplanarPolyg &&
2175                    !coplanarPolyg->isPolyLink( polygon._links[ iL2 ]) &&
2176                    coplanarPolyg < & _polygons[ nbQuadPolygons ])
2177                 break;
2178             if ( iL2 == polygon._links.size() )
2179               coplanarPolyg = 0;
2180           }
2181           if ( 0 /*coplanarPolyg*/ ) // coplanar polygon found
2182           {
2183             freeLinks.resize( freeLinks.size() - polygon._polyLinks.size() );
2184             nbFreeLinks -= polygon._polyLinks.size();
2185
2186             // fill freeLinks with links not shared by coplanarPolyg and polygon
2187             for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
2188               if ( polygon._links[ iL ]._link->_faces[1] &&
2189                    polygon._links[ iL ]._link->_faces[0] != coplanarPolyg )
2190               {
2191                 _Face* p = polygon._links[ iL ]._link->_faces[0];
2192                 for ( size_t iL2 = 0; iL2 < p->_links.size(); ++iL2 )
2193                   if ( p->_links[ iL2 ]._link == polygon._links[ iL ]._link )
2194                   {
2195                     freeLinks.push_back( & p->_links[ iL2 ] );
2196                     ++nbFreeLinks;
2197                     freeLinks.back()->RemoveFace( &polygon );
2198                     break;
2199                   }
2200               }
2201             for ( size_t iL = 0; iL < coplanarPolyg->_links.size(); ++iL )
2202               if ( coplanarPolyg->_links[ iL ]._link->_faces[1] &&
2203                    coplanarPolyg->_links[ iL ]._link->_faces[1] != &polygon )
2204               {
2205                 _Face* p = coplanarPolyg->_links[ iL ]._link->_faces[0];
2206                 if ( p == coplanarPolyg )
2207                   p = coplanarPolyg->_links[ iL ]._link->_faces[1];
2208                 for ( size_t iL2 = 0; iL2 < p->_links.size(); ++iL2 )
2209                   if ( p->_links[ iL2 ]._link == coplanarPolyg->_links[ iL ]._link )
2210                   {
2211                     freeLinks.push_back( & p->_links[ iL2 ] );
2212                     ++nbFreeLinks;
2213                     freeLinks.back()->RemoveFace( coplanarPolyg );
2214                     break;
2215                   }
2216               }
2217             // set coplanarPolyg to be re-created next
2218             for ( size_t iP = 0; iP < _polygons.size(); ++iP )
2219               if ( coplanarPolyg == & _polygons[ iP ] )
2220               {
2221                 iPolygon = iP;
2222                 _polygons[ iPolygon ]._links.clear();
2223                 _polygons[ iPolygon ]._polyLinks.clear();
2224                 break;
2225               }
2226             if ( freeLinks.back() == &polygon._links.back() )
2227             {
2228               freeLinks.pop_back();
2229               --nbFreeLinks;
2230             }
2231             _polygons.pop_back();
2232             usedFaceIDs.erase( curFace );
2233             continue;
2234           } // if ( coplanarPolyg )
2235         } // if ( hasEdgeIntersections )
2236
2237         iPolygon = _polygons.size();
2238
2239       } // end of case ( polygon._links.size() > 2 )
2240     } // while ( nbFreeLinks > 0 )
2241
2242     if ( ! checkPolyhedronSize() )
2243     {
2244       return;
2245     }
2246
2247     // create a classic cell if possible
2248     const int nbNodes = _nbCornerNodes + nbIntersections;
2249     bool isClassicElem = false;
2250     if (      nbNodes == 8 && _polygons.size() == 6 ) isClassicElem = addHexa();
2251     else if ( nbNodes == 4 && _polygons.size() == 4 ) isClassicElem = addTetra();
2252     else if ( nbNodes == 6 && _polygons.size() == 5 ) isClassicElem = addPenta();
2253     else if ( nbNodes == 5 && _polygons.size() == 5 ) isClassicElem = addPyra ();
2254     if ( !isClassicElem )
2255     {
2256       _volumeDefs._nodes.clear();
2257       _volumeDefs._quantities.clear();
2258
2259       for ( size_t iF = 0; iF < _polygons.size(); ++iF )
2260       {
2261         const size_t nbLinks = _polygons[ iF ]._links.size();
2262         _volumeDefs._quantities.push_back( nbLinks );
2263         for ( size_t iL = 0; iL < nbLinks; ++iL )
2264           _volumeDefs._nodes.push_back( _polygons[ iF ]._links[ iL ].FirstNode() );
2265       }
2266     }
2267   }
2268   //================================================================================
2269   /*!
2270    * \brief Create elements in the mesh
2271    */
2272   int Hexahedron::MakeElements(SMESH_MesherHelper&                      helper,
2273                                const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap)
2274   {
2275     SMESHDS_Mesh* mesh = helper.GetMeshDS();
2276
2277     size_t nbCells[3] = { _grid->_coords[0].size() - 1,
2278                           _grid->_coords[1].size() - 1,
2279                           _grid->_coords[2].size() - 1 };
2280     const size_t nbGridCells = nbCells[0] * nbCells[1] * nbCells[2];
2281     vector< Hexahedron* > intersectedHex( nbGridCells, 0 );
2282     int nbIntHex = 0;
2283
2284     // set intersection nodes from GridLine's to links of intersectedHex
2285     int i,j,k, iDirOther[3][2] = {{ 1,2 },{ 0,2 },{ 0,1 }};
2286     for ( int iDir = 0; iDir < 3; ++iDir )
2287     {
2288       int dInd[4][3] = { {0,0,0}, {0,0,0}, {0,0,0}, {0,0,0} };
2289       dInd[1][ iDirOther[iDir][0] ] = -1;
2290       dInd[2][ iDirOther[iDir][1] ] = -1;
2291       dInd[3][ iDirOther[iDir][0] ] = -1; dInd[3][ iDirOther[iDir][1] ] = -1;
2292       // loop on GridLine's parallel to iDir
2293       LineIndexer lineInd = _grid->GetLineIndexer( iDir );
2294       for ( ; lineInd.More(); ++lineInd )
2295       {
2296         GridLine& line = _grid->_lines[ iDir ][ lineInd.LineIndex() ];
2297         multiset< F_IntersectPoint >::const_iterator ip = line._intPoints.begin();
2298         for ( ; ip != line._intPoints.end(); ++ip )
2299         {
2300           // if ( !ip->_node ) continue; // intersection at a grid node
2301           lineInd.SetIndexOnLine( ip->_indexOnLine );
2302           for ( int iL = 0; iL < 4; ++iL ) // loop on 4 cells sharing a link
2303           {
2304             i = int(lineInd.I()) + dInd[iL][0];
2305             j = int(lineInd.J()) + dInd[iL][1];
2306             k = int(lineInd.K()) + dInd[iL][2];
2307             if ( i < 0 || i >= nbCells[0] ||
2308                  j < 0 || j >= nbCells[1] ||
2309                  k < 0 || k >= nbCells[2] ) continue;
2310
2311             const size_t hexIndex = _grid->CellIndex( i,j,k );
2312             Hexahedron *& hex = intersectedHex[ hexIndex ];
2313             if ( !hex)
2314             {
2315               hex = new Hexahedron( *this );
2316               hex->_i = i;
2317               hex->_j = j;
2318               hex->_k = k;
2319               ++nbIntHex;
2320             }
2321             const int iLink = iL + iDir * 4;
2322             hex->_hexLinks[iLink]._fIntPoints.push_back( &(*ip) );
2323             hex->_nbFaceIntNodes += bool( ip->_node );
2324           }
2325         }
2326       }
2327     }
2328
2329     // implement geom edges into the mesh
2330     addEdges( helper, intersectedHex, edge2faceIDsMap );
2331
2332     // add not split hexadrons to the mesh
2333     int nbAdded = 0;
2334     vector<int> intHexInd( nbIntHex );
2335     nbIntHex = 0;
2336     for ( size_t i = 0; i < intersectedHex.size(); ++i )
2337     {
2338       Hexahedron * & hex = intersectedHex[ i ];
2339       if ( hex )
2340       {
2341         intHexInd[ nbIntHex++ ] = i;
2342         if ( hex->_nbFaceIntNodes > 0 || hex->_eIntPoints.size() > 0 )
2343           continue; // treat intersected hex later
2344         this->init( hex->_i, hex->_j, hex->_k );
2345       }
2346       else
2347       {    
2348         this->init( i );
2349       }
2350       if (( _nbCornerNodes == 8 ) &&
2351           ( _nbBndNodes < _nbCornerNodes || !isInHole() ))
2352       {
2353         // order of _hexNodes is defined by enum SMESH_Block::TShapeID
2354         SMDS_MeshElement* el =
2355           mesh->AddVolume( _hexNodes[0].Node(), _hexNodes[2].Node(),
2356                            _hexNodes[3].Node(), _hexNodes[1].Node(),
2357                            _hexNodes[4].Node(), _hexNodes[6].Node(),
2358                            _hexNodes[7].Node(), _hexNodes[5].Node() );
2359         mesh->SetMeshElementOnShape( el, helper.GetSubShapeID() );
2360         ++nbAdded;
2361         if ( hex )
2362         {
2363           delete hex;
2364           intersectedHex[ i ] = 0;
2365           --nbIntHex;
2366         }
2367       }
2368       else if ( _nbCornerNodes > 3  && !hex )
2369       {
2370         // all intersection of hex with geometry are at grid nodes
2371         hex = new Hexahedron( *this );
2372         hex->_i = _i;
2373         hex->_j = _j;
2374         hex->_k = _k;
2375         intHexInd.push_back(0);
2376         intHexInd[ nbIntHex++ ] = i;
2377       }
2378     }
2379
2380     // add elements resulted from hexadron intersection
2381 #ifdef WITH_TBB
2382     intHexInd.resize( nbIntHex );
2383     tbb::parallel_for ( tbb::blocked_range<size_t>( 0, nbIntHex ),
2384                         ParallelHexahedron( intersectedHex, intHexInd ),
2385                         tbb::simple_partitioner()); // ComputeElements() is called here
2386     for ( size_t i = 0; i < intHexInd.size(); ++i )
2387       if ( Hexahedron * hex = intersectedHex[ intHexInd[ i ]] )
2388         nbAdded += hex->addElements( helper );
2389 #else
2390     for ( size_t i = 0; i < intHexInd.size(); ++i )
2391       if ( Hexahedron * hex = intersectedHex[ intHexInd[ i ]] )
2392       {
2393         hex->ComputeElements();
2394         nbAdded += hex->addElements( helper );
2395       }
2396 #endif
2397
2398     for ( size_t i = 0; i < intersectedHex.size(); ++i )
2399       if ( intersectedHex[ i ] )
2400         delete intersectedHex[ i ];
2401
2402     return nbAdded;
2403   }
2404
2405   //================================================================================
2406   /*!
2407    * \brief Implements geom edges into the mesh
2408    */
2409   void Hexahedron::addEdges(SMESH_MesherHelper&                      helper,
2410                             vector< Hexahedron* >&                   hexes,
2411                             const map< TGeomID, vector< TGeomID > >& edge2faceIDsMap)
2412   {
2413     if ( edge2faceIDsMap.empty() ) return;
2414
2415     // Prepare planes for intersecting with EDGEs
2416     GridPlanes pln[3];
2417     {
2418       for ( int iDirZ = 0; iDirZ < 3; ++iDirZ ) // iDirZ gives normal direction to planes
2419       {
2420         GridPlanes& planes = pln[ iDirZ ];
2421         int iDirX = ( iDirZ + 1 ) % 3;
2422         int iDirY = ( iDirZ + 2 ) % 3;
2423         planes._zNorm  = ( _grid->_axes[ iDirX ] ^ _grid->_axes[ iDirY ] ).Normalized();
2424         planes._zProjs.resize ( _grid->_coords[ iDirZ ].size() );
2425         planes._zProjs [0] = 0;
2426         const double       zFactor = _grid->_axes[ iDirZ ] * planes._zNorm;
2427         const vector< double > & u = _grid->_coords[ iDirZ ];
2428         for ( int i = 1; i < planes._zProjs.size(); ++i )
2429         {
2430           planes._zProjs [i] = zFactor * ( u[i] - u[0] );
2431         }
2432       }
2433     }
2434     const double deflection = _grid->_minCellSize / 20.;
2435     const double tol        = _grid->_tol;
2436     E_IntersectPoint ip;
2437
2438     // Intersect EDGEs with the planes
2439     map< TGeomID, vector< TGeomID > >::const_iterator e2fIt = edge2faceIDsMap.begin();
2440     for ( ; e2fIt != edge2faceIDsMap.end(); ++e2fIt )
2441     {
2442       const TGeomID  edgeID = e2fIt->first;
2443       const TopoDS_Edge & E = TopoDS::Edge( _grid->_shapes( edgeID ));
2444       BRepAdaptor_Curve curve( E );
2445       TopoDS_Vertex v1 = helper.IthVertex( 0, E, false ); 
2446       TopoDS_Vertex v2 = helper.IthVertex( 1, E, false ); 
2447
2448       ip._faceIDs = e2fIt->second;
2449       ip._shapeID = edgeID;
2450
2451       // discretize the EGDE
2452       GCPnts_UniformDeflection discret( curve, deflection, true );
2453       if ( !discret.IsDone() || discret.NbPoints() < 2 )
2454         continue;
2455
2456       // perform intersection
2457       for ( int iDirZ = 0; iDirZ < 3; ++iDirZ )
2458       {
2459         GridPlanes& planes = pln[ iDirZ ];
2460         int      iDirX = ( iDirZ + 1 ) % 3;
2461         int      iDirY = ( iDirZ + 2 ) % 3;
2462         double    xLen = _grid->_coords[ iDirX ].back() - _grid->_coords[ iDirX ][0];
2463         double    yLen = _grid->_coords[ iDirY ].back() - _grid->_coords[ iDirY ][0];
2464         double    zLen = _grid->_coords[ iDirZ ].back() - _grid->_coords[ iDirZ ][0];
2465         int dIJK[3], d000[3] = { 0,0,0 };
2466         double o[3] = { _grid->_coords[0][0],
2467                         _grid->_coords[1][0],
2468                         _grid->_coords[2][0] };
2469
2470         // locate the 1st point of a segment within the grid
2471         gp_XYZ p1     = discret.Value( 1 ).XYZ();
2472         double u1     = discret.Parameter( 1 );
2473         double zProj1 = planes._zNorm * ( p1 - _grid->_origin );
2474
2475         _grid->ComputeUVW( p1, ip._uvw );
2476         int iX1 = int(( ip._uvw[iDirX] - o[iDirX]) / xLen * (_grid->_coords[ iDirX ].size() - 1));
2477         int iY1 = int(( ip._uvw[iDirY] - o[iDirY]) / yLen * (_grid->_coords[ iDirY ].size() - 1));
2478         int iZ1 = int(( ip._uvw[iDirZ] - o[iDirZ]) / zLen * (_grid->_coords[ iDirZ ].size() - 1));
2479         locateValue( iX1, ip._uvw[iDirX], _grid->_coords[ iDirX ], dIJK[ iDirX ], tol );
2480         locateValue( iY1, ip._uvw[iDirY], _grid->_coords[ iDirY ], dIJK[ iDirY ], tol );
2481         locateValue( iZ1, ip._uvw[iDirZ], _grid->_coords[ iDirZ ], dIJK[ iDirZ ], tol );
2482
2483         int ijk[3]; // grid index where a segment intersect a plane
2484         ijk[ iDirX ] = iX1;
2485         ijk[ iDirY ] = iY1;
2486         ijk[ iDirZ ] = iZ1;
2487
2488         // add the 1st vertex point to a hexahedron
2489         if ( iDirZ == 0 )
2490         {
2491           ip._point   = p1;
2492           ip._shapeID = _grid->_shapes.Add( v1 );
2493           _grid->_edgeIntP.push_back( ip );
2494           if ( !addIntersection( _grid->_edgeIntP.back(), hexes, ijk, d000 ))
2495             _grid->_edgeIntP.pop_back();
2496           ip._shapeID = edgeID;
2497         }
2498         for ( int iP = 2; iP <= discret.NbPoints(); ++iP )
2499         {
2500           // locate the 2nd point of a segment within the grid
2501           gp_XYZ p2     = discret.Value( iP ).XYZ();
2502           double u2     = discret.Parameter( iP );
2503           double zProj2 = planes._zNorm * ( p2 - _grid->_origin );
2504           int    iZ2    = iZ1;
2505           if ( Abs( zProj2 - zProj1 ) <= std::numeric_limits<double>::min() )
2506             continue;
2507           locateValue( iZ2, zProj2, planes._zProjs, dIJK[ iDirZ ], tol );
2508
2509           // treat intersections with planes between 2 end points of a segment
2510           int dZ = ( iZ1 <= iZ2 ) ? +1 : -1;
2511           int iZ = iZ1 + ( iZ1 < iZ2 );
2512           for ( int i = 0, nb = Abs( iZ1 - iZ2 ); i < nb; ++i, iZ += dZ )
2513           {
2514             ip._point = findIntPoint( u1, zProj1, u2, zProj2,
2515                                       planes._zProjs[ iZ ],
2516                                       curve, planes._zNorm, _grid->_origin );
2517             _grid->ComputeUVW( ip._point.XYZ(), ip._uvw );
2518             locateValue( ijk[iDirX], ip._uvw[iDirX], _grid->_coords[iDirX], dIJK[iDirX], tol );
2519             locateValue( ijk[iDirY], ip._uvw[iDirY], _grid->_coords[iDirY], dIJK[iDirY], tol );
2520             ijk[ iDirZ ] = iZ;
2521
2522             // add ip to hex "above" the plane
2523             _grid->_edgeIntP.push_back( ip );
2524             dIJK[ iDirZ ] = 0;
2525             bool added = addIntersection(_grid->_edgeIntP.back(), hexes, ijk, dIJK);
2526
2527             // add ip to hex "below" the plane
2528             ijk[ iDirZ ] = iZ-1;
2529             if ( !addIntersection( _grid->_edgeIntP.back(), hexes, ijk, dIJK ) &&
2530                  !added)
2531               _grid->_edgeIntP.pop_back();
2532           }
2533           iZ1    = iZ2;
2534           p1     = p2;
2535           u1     = u2;
2536           zProj1 = zProj2;
2537         }
2538         // add the 2nd vertex point to a hexahedron
2539         if ( iDirZ == 0 )
2540         {
2541           ip._shapeID = _grid->_shapes.Add( v2 );
2542           ip._point = p1;
2543           _grid->ComputeUVW( p1, ip._uvw );
2544           locateValue( ijk[iDirX], ip._uvw[iDirX], _grid->_coords[iDirX], dIJK[iDirX], tol );
2545           locateValue( ijk[iDirY], ip._uvw[iDirY], _grid->_coords[iDirY], dIJK[iDirY], tol );
2546           ijk[ iDirZ ] = iZ1;
2547           _grid->_edgeIntP.push_back( ip );
2548           if ( !addIntersection( _grid->_edgeIntP.back(), hexes, ijk, d000 ))
2549             _grid->_edgeIntP.pop_back();
2550           ip._shapeID = edgeID;
2551         }
2552       } // loop on 3 grid directions
2553     } // loop on EDGEs
2554
2555   }
2556
2557   //================================================================================
2558   /*!
2559    * \brief Finds intersection of a curve with a plane
2560    *  \param [in] u1 - parameter of one curve point
2561    *  \param [in] proj1 - projection of the curve point to the plane normal
2562    *  \param [in] u2 - parameter of another curve point
2563    *  \param [in] proj2 - projection of the other curve point to the plane normal
2564    *  \param [in] proj - projection of a point where the curve intersects the plane
2565    *  \param [in] curve - the curve
2566    *  \param [in] axis - the plane normal
2567    *  \param [in] origin - the plane origin
2568    *  \return gp_Pnt - the found intersection point
2569    */
2570   gp_Pnt Hexahedron::findIntPoint( double u1, double proj1,
2571                                    double u2, double proj2,
2572                                    double proj,
2573                                    BRepAdaptor_Curve& curve,
2574                                    const gp_XYZ& axis,
2575                                    const gp_XYZ& origin)
2576   {
2577     double r = (( proj - proj1 ) / ( proj2 - proj1 ));
2578     double u = u1 * ( 1 - r ) + u2 * r;
2579     gp_Pnt p = curve.Value( u );
2580     double newProj =  axis * ( p.XYZ() - origin );
2581     if ( Abs( proj - newProj ) > _grid->_tol / 10. )
2582     {
2583       if ( r > 0.5 )
2584         return findIntPoint( u2, proj2, u, newProj, proj, curve, axis, origin );
2585       else
2586         return findIntPoint( u1, proj2, u, newProj, proj, curve, axis, origin );
2587     }
2588     return p;
2589   }
2590
2591   //================================================================================
2592   /*!
2593    * \brief Returns indices of a hexahedron sub-entities holding a point
2594    *  \param [in] ip - intersection point
2595    *  \param [out] facets - 0-3 facets holding a point
2596    *  \param [out] sub - index of a vertex or an edge holding a point
2597    *  \return int - number of facets holding a point
2598    */
2599   int Hexahedron::getEntity( const E_IntersectPoint* ip, int* facets, int& sub )
2600   {
2601     enum { X = 1, Y = 2, Z = 4 }; // == 001, 010, 100
2602     int nbFacets = 0;
2603     int vertex = 0, egdeMask = 0;
2604
2605     if ( Abs( _grid->_coords[0][ _i   ] - ip->_uvw[0] ) < _grid->_tol ) {
2606       facets[ nbFacets++ ] = SMESH_Block::ID_F0yz;
2607       egdeMask |= X;
2608     }
2609     else if ( Abs( _grid->_coords[0][ _i+1 ] - ip->_uvw[0] ) < _grid->_tol ) {
2610       facets[ nbFacets++ ] = SMESH_Block::ID_F1yz;
2611       vertex   |= X;
2612       egdeMask |= X;
2613     }
2614     if ( Abs( _grid->_coords[1][ _j   ] - ip->_uvw[1] ) < _grid->_tol ) {
2615       facets[ nbFacets++ ] = SMESH_Block::ID_Fx0z;
2616       egdeMask |= Y;
2617     }
2618     else if ( Abs( _grid->_coords[1][ _j+1 ] - ip->_uvw[1] ) < _grid->_tol ) {
2619       facets[ nbFacets++ ] = SMESH_Block::ID_Fx1z;
2620       vertex   |= Y;
2621       egdeMask |= Y;
2622     }
2623     if ( Abs( _grid->_coords[2][ _k   ] - ip->_uvw[2] ) < _grid->_tol ) {
2624       facets[ nbFacets++ ] = SMESH_Block::ID_Fxy0;
2625       egdeMask |= Z;
2626     }
2627     else if ( Abs( _grid->_coords[2][ _k+1 ] - ip->_uvw[2] ) < _grid->_tol ) {
2628       facets[ nbFacets++ ] = SMESH_Block::ID_Fxy1;
2629       vertex   |= Z;
2630       egdeMask |= Z;
2631     }
2632
2633     switch ( nbFacets )
2634     {
2635     case 0: sub = 0;         break;
2636     case 1: sub = facets[0]; break;
2637     case 2: {
2638       const int edge [3][8] = {
2639         { SMESH_Block::ID_E00z, SMESH_Block::ID_E10z,
2640           SMESH_Block::ID_E01z, SMESH_Block::ID_E11z },
2641         { SMESH_Block::ID_E0y0, SMESH_Block::ID_E1y0, 0, 0,
2642           SMESH_Block::ID_E0y1, SMESH_Block::ID_E1y1 },
2643         { SMESH_Block::ID_Ex00, 0, SMESH_Block::ID_Ex10, 0,
2644           SMESH_Block::ID_Ex01, 0, SMESH_Block::ID_Ex11 }
2645       };
2646       switch ( egdeMask ) {
2647       case X | Y: sub = edge[ 0 ][ vertex ]; break;
2648       case X | Z: sub = edge[ 1 ][ vertex ]; break;
2649       default:    sub = edge[ 2 ][ vertex ];
2650       }
2651       break;
2652     }
2653     //case 3:
2654     default:
2655       sub = vertex + SMESH_Block::ID_FirstV;
2656     }
2657
2658     return nbFacets;
2659   }
2660   //================================================================================
2661   /*!
2662    * \brief Adds intersection with an EDGE
2663    */
2664   bool Hexahedron::addIntersection( const E_IntersectPoint& ip,
2665                                     vector< Hexahedron* >&  hexes,
2666                                     int ijk[], int dIJK[] )
2667   {
2668     bool added = false;
2669
2670     size_t hexIndex[4] = {
2671       _grid->CellIndex( ijk[0], ijk[1], ijk[2] ),
2672       dIJK[0] ? _grid->CellIndex( ijk[0]+dIJK[0], ijk[1], ijk[2] ) : -1,
2673       dIJK[1] ? _grid->CellIndex( ijk[0], ijk[1]+dIJK[1], ijk[2] ) : -1,
2674       dIJK[2] ? _grid->CellIndex( ijk[0], ijk[1], ijk[2]+dIJK[2] ) : -1
2675     };
2676     for ( int i = 0; i < 4; ++i )
2677     {
2678       if ( /*0 <= hexIndex[i] &&*/ hexIndex[i] < hexes.size() && hexes[ hexIndex[i] ] )
2679       {
2680         Hexahedron* h = hexes[ hexIndex[i] ];
2681         // check if ip is really inside the hex
2682 #ifdef _DEBUG_
2683         if (( _grid->_coords[0][ h->_i   ] - _grid->_tol > ip._uvw[0] ) ||
2684             ( _grid->_coords[0][ h->_i+1 ] + _grid->_tol < ip._uvw[0] ) ||
2685             ( _grid->_coords[1][ h->_j   ] - _grid->_tol > ip._uvw[1] ) ||
2686             ( _grid->_coords[1][ h->_j+1 ] + _grid->_tol < ip._uvw[1] ) ||
2687             ( _grid->_coords[2][ h->_k   ] - _grid->_tol > ip._uvw[2] ) ||
2688             ( _grid->_coords[2][ h->_k+1 ] + _grid->_tol < ip._uvw[2] ))
2689           throw SALOME_Exception("ip outside a hex");
2690 #endif
2691         h->_eIntPoints.push_back( & ip );
2692         added = true;
2693       }
2694     }
2695     return added;
2696   }
2697   //================================================================================
2698   /*!
2699    * \brief Finds nodes at a path from one node to another via intersections with EDGEs
2700    */
2701   bool Hexahedron::findChain( _Node*          n1,
2702                               _Node*          n2,
2703                               _Face&          quad,
2704                               vector<_Node*>& chn )
2705   {
2706     chn.clear();
2707     chn.push_back( n1 );
2708     for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
2709       if ( !quad._eIntNodes[ iP ]->IsUsedInFace( &quad ) &&
2710            n1->IsLinked( quad._eIntNodes[ iP ]->_intPoint ) &&
2711            n2->IsLinked( quad._eIntNodes[ iP ]->_intPoint ))
2712       {
2713         chn.push_back( quad._eIntNodes[ iP ]);
2714         chn.push_back( n2 );
2715         quad._eIntNodes[ iP ]->_usedInFace = &quad;
2716         return true;
2717       }
2718     bool found;
2719     do
2720     {
2721       found = false;
2722       for ( size_t iP = 0; iP < quad._eIntNodes.size(); ++iP )
2723         if ( !quad._eIntNodes[ iP ]->IsUsedInFace( &quad ) &&
2724              chn.back()->IsLinked( quad._eIntNodes[ iP ]->_intPoint ))
2725         {
2726           chn.push_back( quad._eIntNodes[ iP ]);
2727           found = quad._eIntNodes[ iP ]->_usedInFace = &quad;
2728           break;
2729         }
2730     } while ( found && ! chn.back()->IsLinked( n2->_intPoint ) );
2731
2732     if ( chn.back() != n2 && chn.back()->IsLinked( n2->_intPoint ))
2733       chn.push_back( n2 );
2734
2735     return chn.size() > 1;
2736   }
2737   //================================================================================
2738   /*!
2739    * \brief Try to heal a polygon whose ends are not connected
2740    */
2741   bool Hexahedron::closePolygon( _Face* polygon, vector<_Node*>& chainNodes ) const
2742   {
2743     int i = -1, nbLinks = polygon->_links.size();
2744     if ( nbLinks < 3 )
2745       return false;
2746     vector< _OrientedLink > newLinks;
2747     // find a node lying on the same FACE as the last one
2748     _Node*   node = polygon->_links.back().LastNode();
2749     int avoidFace = node->IsLinked( polygon->_links.back().FirstNode()->_intPoint );
2750     for ( i = nbLinks - 2; i >= 0; --i )
2751       if ( node->IsLinked( polygon->_links[i].FirstNode()->_intPoint, avoidFace ))
2752         break;
2753     if ( i >= 0 )
2754     {
2755       for ( ; i < nbLinks; ++i )
2756         newLinks.push_back( polygon->_links[i] );
2757     }
2758     else
2759     {
2760       // find a node lying on the same FACE as the first one
2761       node      = polygon->_links[0].FirstNode();
2762       avoidFace = node->IsLinked( polygon->_links[0].LastNode()->_intPoint );
2763       for ( i = 1; i < nbLinks; ++i )
2764         if ( node->IsLinked( polygon->_links[i].LastNode()->_intPoint, avoidFace ))
2765           break;
2766       if ( i < nbLinks )
2767         for ( nbLinks = i + 1, i = 0; i < nbLinks; ++i )
2768           newLinks.push_back( polygon->_links[i] );
2769     }
2770     if ( newLinks.size() > 1 )
2771     {
2772       polygon->_links.swap( newLinks );
2773       chainNodes.clear();
2774       chainNodes.push_back( polygon->_links.back().LastNode() );
2775       chainNodes.push_back( polygon->_links[0].FirstNode() );
2776       return true;
2777     }
2778     return false;
2779   }
2780   //================================================================================
2781   /*!
2782    * \brief Checks transition at the 1st node of a link
2783    */
2784   bool Hexahedron::is1stNodeOut( _Link& link /*int iLink*/ ) const
2785   {
2786     // new version is for the case: tangent transition at the 1st node
2787     bool isOut = false;
2788     if ( link._fIntNodes.size() > 1 )
2789     {
2790       // check transition at the next intersection
2791       switch ( link._fIntPoints[1]->_transition ) {
2792       case Trans_OUT: return false;
2793       case Trans_IN : return true;
2794       default: ; // tangent transition
2795       }
2796     }
2797     gp_Pnt p1 = link._nodes[0]->Point();
2798     gp_Pnt p2 = link._nodes[1]->Point();
2799     gp_Pnt testPnt = 0.8 * p1.XYZ() + 0.2 * p2.XYZ();
2800
2801     TGeomID          faceID = link._fIntPoints[0]->_faceIDs[0];
2802     const TopoDS_Face& face = TopoDS::Face( _grid->_shapes( faceID ));
2803     TopLoc_Location loc;
2804     GeomAPI_ProjectPointOnSurf& proj =
2805       _grid->_helper->GetProjector( face, loc, 0.1*_grid->_tol );
2806     testPnt.Transform( loc );
2807     proj.Perform( testPnt );
2808     if ( proj.IsDone() &&
2809          proj.NbPoints() > 0 &&
2810          proj.LowerDistance() > _grid->_tol )
2811     {
2812       Quantity_Parameter u,v;
2813       proj.LowerDistanceParameters( u,v );
2814       gp_Dir normal;
2815       if ( GeomLib::NormEstim( BRep_Tool::Surface( face, loc ),
2816                                gp_Pnt2d( u,v ),
2817                                0.1*_grid->_tol,
2818                                normal ) < 3 )
2819       {
2820         if ( face.Orientation() == TopAbs_REVERSED )
2821           normal.Reverse();
2822         gp_Vec v( proj.NearestPoint(), testPnt );
2823         return v * normal > 0;
2824       }
2825     }
2826     return isOut;
2827   }
2828   //================================================================================
2829   /*!
2830    * \brief Adds computed elements to the mesh
2831    */
2832   int Hexahedron::addElements(SMESH_MesherHelper& helper)
2833   {
2834     int nbAdded = 0;
2835     // add elements resulted from hexahedron intersection
2836     //for ( size_t i = 0; i < _volumeDefs.size(); ++i )
2837     {
2838       vector< const SMDS_MeshNode* > nodes( _volumeDefs._nodes.size() );
2839       for ( size_t iN = 0; iN < nodes.size(); ++iN )
2840         if ( !( nodes[iN] = _volumeDefs._nodes[iN]->Node() ))
2841         {
2842           if ( const E_IntersectPoint* eip = _volumeDefs._nodes[iN]->EdgeIntPnt() )
2843             nodes[iN] = _volumeDefs._nodes[iN]->_intPoint->_node =
2844               helper.AddNode( eip->_point.X(),
2845                               eip->_point.Y(),
2846                               eip->_point.Z() );
2847           else
2848             throw SALOME_Exception("Bug: no node at intersection point");
2849         }
2850
2851       if ( !_volumeDefs._quantities.empty() )
2852       {
2853         helper.AddPolyhedralVolume( nodes, _volumeDefs._quantities );
2854       }
2855       else
2856       {
2857         switch ( nodes.size() )
2858         {
2859         case 8: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],
2860                                   nodes[4],nodes[5],nodes[6],nodes[7] );
2861           break;
2862         case 4: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3] );
2863           break;
2864         case 6: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3], nodes[4],nodes[5] );
2865           break;
2866         case 5:
2867           helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],nodes[4] );
2868           break;
2869         }
2870       }
2871       nbAdded += int ( _volumeDefs._nodes.size() > 0 );
2872     }
2873
2874     return nbAdded;
2875   }
2876   //================================================================================
2877   /*!
2878    * \brief Return true if the element is in a hole
2879    */
2880   bool Hexahedron::isInHole() const
2881   {
2882     if ( !_vIntNodes.empty() )
2883       return false;
2884
2885     const int ijk[3] = { _i, _j, _k };
2886     F_IntersectPoint curIntPnt;
2887
2888     // consider a cell to be in a hole if all links in any direction
2889     // comes OUT of geometry
2890     for ( int iDir = 0; iDir < 3; ++iDir )
2891     {
2892       const vector<double>& coords = _grid->_coords[ iDir ];
2893       LineIndexer               li = _grid->GetLineIndexer( iDir );
2894       li.SetIJK( _i,_j,_k );
2895       size_t lineIndex[4] = { li.LineIndex  (),
2896                               li.LineIndex10(),
2897                               li.LineIndex01(),
2898                               li.LineIndex11() };
2899       bool allLinksOut = true, hasLinks = false;
2900       for ( int iL = 0; iL < 4 && allLinksOut; ++iL ) // loop on 4 links parallel to iDir
2901       {
2902         const _Link& link = _hexLinks[ iL + 4*iDir ];
2903         // check transition of the first node of a link
2904         const F_IntersectPoint* firstIntPnt = 0;
2905         if ( link._nodes[0]->Node() ) // 1st node is a hexa corner
2906         {
2907           curIntPnt._paramOnLine = coords[ ijk[ iDir ]] - coords[0];
2908           const GridLine& line = _grid->_lines[ iDir ][ lineIndex[ iL ]];
2909           multiset< F_IntersectPoint >::const_iterator ip =
2910             line._intPoints.upper_bound( curIntPnt );
2911           --ip;
2912           firstIntPnt = &(*ip);
2913         }
2914         else if ( !link._fIntPoints.empty() )
2915         {
2916           firstIntPnt = link._fIntPoints[0];
2917         }
2918
2919         if ( firstIntPnt )
2920         {
2921           hasLinks = true;
2922           allLinksOut = ( firstIntPnt->_transition == Trans_OUT );
2923         }
2924       }
2925       if ( hasLinks && allLinksOut )
2926         return true;
2927     }
2928     return false;
2929   }
2930
2931   //================================================================================
2932   /*!
2933    * \brief Return true if a polyhedron passes _sizeThreshold criterion
2934    */
2935   bool Hexahedron::checkPolyhedronSize() const
2936   {
2937     double volume = 0;
2938     for ( size_t iP = 0; iP < _polygons.size(); ++iP )
2939     {
2940       const _Face& polygon = _polygons[iP];
2941       gp_XYZ area (0,0,0);
2942       gp_XYZ p1 = polygon._links[ 0 ].FirstNode()->Point().XYZ();
2943       for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
2944       {
2945         gp_XYZ p2 = polygon._links[ iL ].LastNode()->Point().XYZ();
2946         area += p1 ^ p2;
2947         p1 = p2;
2948       }
2949       volume += p1 * area;
2950     }
2951     volume /= 6;
2952
2953     double initVolume = _sideLength[0] * _sideLength[1] * _sideLength[2];
2954
2955     return volume > initVolume / _sizeThreshold;
2956   }
2957   //================================================================================
2958   /*!
2959    * \brief Tries to create a hexahedron
2960    */
2961   bool Hexahedron::addHexa()
2962   {
2963     if ( _polygons[0]._links.size() != 4 ||
2964          _polygons[1]._links.size() != 4 ||
2965          _polygons[2]._links.size() != 4 ||
2966          _polygons[3]._links.size() != 4 ||
2967          _polygons[4]._links.size() != 4 ||
2968          _polygons[5]._links.size() != 4   )
2969       return false;
2970     _Node* nodes[8];
2971     int nbN = 0;
2972     for ( int iL = 0; iL < 4; ++iL )
2973     {
2974       // a base node
2975       nodes[iL] = _polygons[0]._links[iL].FirstNode();
2976       ++nbN;
2977
2978       // find a top node above the base node
2979       _Link* link = _polygons[0]._links[iL]._link;
2980       if ( !link->_faces[0] || !link->_faces[1] )
2981         return debugDumpLink( link );
2982       // a quadrangle sharing <link> with _polygons[0]
2983       _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[0] )];
2984       for ( int i = 0; i < 4; ++i )
2985         if ( quad->_links[i]._link == link )
2986         {
2987           // 1st node of a link opposite to <link> in <quad>
2988           nodes[iL+4] = quad->_links[(i+2)%4].FirstNode();
2989           ++nbN;
2990           break;
2991         }
2992     }
2993     if ( nbN == 8 )
2994       _volumeDefs.set( vector< _Node* >( nodes, nodes+8 ));
2995
2996     return nbN == 8;
2997   }
2998   //================================================================================
2999   /*!
3000    * \brief Tries to create a tetrahedron
3001    */
3002   bool Hexahedron::addTetra()
3003   {
3004     _Node* nodes[4];
3005     nodes[0] = _polygons[0]._links[0].FirstNode();
3006     nodes[1] = _polygons[0]._links[1].FirstNode();
3007     nodes[2] = _polygons[0]._links[2].FirstNode();
3008
3009     _Link* link = _polygons[0]._links[0]._link;
3010     if ( !link->_faces[0] || !link->_faces[1] )
3011       return debugDumpLink( link );
3012
3013     // a triangle sharing <link> with _polygons[0]
3014     _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[0] )];
3015     for ( int i = 0; i < 3; ++i )
3016       if ( tria->_links[i]._link == link )
3017       {
3018         nodes[3] = tria->_links[(i+1)%3].LastNode();
3019         _volumeDefs.set( vector< _Node* >( nodes, nodes+4 ));
3020         return true;
3021       }
3022
3023     return false;
3024   }
3025   //================================================================================
3026   /*!
3027    * \brief Tries to create a pentahedron
3028    */
3029   bool Hexahedron::addPenta()
3030   {
3031     // find a base triangular face
3032     int iTri = -1;
3033     for ( int iF = 0; iF < 5 && iTri < 0; ++iF )
3034       if ( _polygons[ iF ]._links.size() == 3 )
3035         iTri = iF;
3036     if ( iTri < 0 ) return false;
3037
3038     // find nodes
3039     _Node* nodes[6];
3040     int nbN = 0;
3041     for ( int iL = 0; iL < 3; ++iL )
3042     {
3043       // a base node
3044       nodes[iL] = _polygons[ iTri ]._links[iL].FirstNode();
3045       ++nbN;
3046
3047       // find a top node above the base node
3048       _Link* link = _polygons[ iTri ]._links[iL]._link;
3049       if ( !link->_faces[0] || !link->_faces[1] )
3050         return debugDumpLink( link );
3051       // a quadrangle sharing <link> with a base triangle
3052       _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[ iTri ] )];
3053       if ( quad->_links.size() != 4 ) return false;
3054       for ( int i = 0; i < 4; ++i )
3055         if ( quad->_links[i]._link == link )
3056         {
3057           // 1st node of a link opposite to <link> in <quad>
3058           nodes[iL+3] = quad->_links[(i+2)%4].FirstNode();
3059           ++nbN;
3060           break;
3061         }
3062     }
3063     if ( nbN == 6 )
3064       _volumeDefs.set( vector< _Node* >( nodes, nodes+6 ));
3065
3066     return ( nbN == 6 );
3067   }
3068   //================================================================================
3069   /*!
3070    * \brief Tries to create a pyramid
3071    */
3072   bool Hexahedron::addPyra()
3073   {
3074     // find a base quadrangle
3075     int iQuad = -1;
3076     for ( int iF = 0; iF < 5 && iQuad < 0; ++iF )
3077       if ( _polygons[ iF ]._links.size() == 4 )
3078         iQuad = iF;
3079     if ( iQuad < 0 ) return false;
3080
3081     // find nodes
3082     _Node* nodes[5];
3083     nodes[0] = _polygons[iQuad]._links[0].FirstNode();
3084     nodes[1] = _polygons[iQuad]._links[1].FirstNode();
3085     nodes[2] = _polygons[iQuad]._links[2].FirstNode();
3086     nodes[3] = _polygons[iQuad]._links[3].FirstNode();
3087
3088     _Link* link = _polygons[iQuad]._links[0]._link;
3089     if ( !link->_faces[0] || !link->_faces[1] )
3090       return debugDumpLink( link );
3091
3092     // a triangle sharing <link> with a base quadrangle
3093     _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[ iQuad ] )];
3094     if ( tria->_links.size() != 3 ) return false;
3095     for ( int i = 0; i < 3; ++i )
3096       if ( tria->_links[i]._link == link )
3097       {
3098         nodes[4] = tria->_links[(i+1)%3].LastNode();
3099         _volumeDefs.set( vector< _Node* >( nodes, nodes+5 ));
3100         return true;
3101       }
3102
3103     return false;
3104   }
3105   //================================================================================
3106   /*!
3107    * \brief Dump a link and return \c false
3108    */
3109   bool Hexahedron::debugDumpLink( Hexahedron::_Link* link )
3110   {
3111 #ifdef _DEBUG_
3112     gp_Pnt p1 = link->_nodes[0]->Point(), p2 = link->_nodes[1]->Point();
3113     cout << "BUG: not shared link. IKJ = ( "<< _i << " " << _j << " " << _k << " )" << endl
3114          << "n1 (" << p1.X() << ", "<< p1.Y() << ", "<< p1.Z() << " )" << endl
3115          << "n2 (" << p2.X() << ", "<< p2.Y() << ", "<< p2.Z() << " )" << endl;
3116 #endif
3117     return false;
3118   }
3119
3120   //================================================================================
3121   /*!
3122    * \brief computes exact bounding box with axes parallel to given ones
3123    */
3124   //================================================================================
3125
3126   void getExactBndBox( const vector< TopoDS_Shape >& faceVec,
3127                        const double*                 axesDirs,
3128                        Bnd_Box&                      shapeBox )
3129   {
3130     BRep_Builder b;
3131     TopoDS_Compound allFacesComp;
3132     b.MakeCompound( allFacesComp );
3133     for ( size_t iF = 0; iF < faceVec.size(); ++iF )
3134       b.Add( allFacesComp, faceVec[ iF ] );
3135
3136     double sP[6]; // aXmin, aYmin, aZmin, aXmax, aYmax, aZmax
3137     shapeBox.Get(sP[0],sP[1],sP[2],sP[3],sP[4],sP[5]);
3138     double farDist = 0;
3139     for ( int i = 0; i < 6; ++i )
3140       farDist = Max( farDist, 10 * sP[i] );
3141
3142     gp_XYZ axis[3] = { gp_XYZ( axesDirs[0], axesDirs[1], axesDirs[2] ),
3143                        gp_XYZ( axesDirs[3], axesDirs[4], axesDirs[5] ),
3144                        gp_XYZ( axesDirs[6], axesDirs[7], axesDirs[8] ) };
3145     axis[0].Normalize();
3146     axis[1].Normalize();
3147     axis[2].Normalize();
3148
3149     gp_Mat basis( axis[0], axis[1], axis[2] );
3150     gp_Mat bi = basis.Inverted();
3151
3152     gp_Pnt pMin, pMax;
3153     for ( int iDir = 0; iDir < 3; ++iDir )
3154     {
3155       gp_XYZ axis0 = axis[ iDir ];
3156       gp_XYZ axis1 = axis[ ( iDir + 1 ) % 3 ];
3157       gp_XYZ axis2 = axis[ ( iDir + 2 ) % 3 ];
3158       for ( int isMax = 0; isMax < 2; ++isMax )
3159       {
3160         double shift = isMax ? farDist : -farDist;
3161         gp_XYZ orig = shift * axis0;
3162         gp_XYZ norm = axis1 ^ axis2;
3163         gp_Pln pln( orig, norm );
3164         norm = pln.Axis().Direction().XYZ();
3165         BRepBuilderAPI_MakeFace plane( pln, -farDist, farDist, -farDist, farDist );
3166
3167         gp_Pnt& pAxis = isMax ? pMax : pMin;
3168         gp_Pnt pPlane, pFaces;
3169         double dist = GEOMUtils::GetMinDistance( plane, allFacesComp, pPlane, pFaces );
3170         if ( dist < 0 )
3171         {
3172           Bnd_B3d bb;
3173           gp_XYZ corner;
3174           for ( int i = 0; i < 2; ++i ) {
3175             corner.SetCoord( 1, sP[ i*3 ]);
3176             for ( int j = 0; j < 2; ++j ) {
3177               corner.SetCoord( 2, sP[ i*3 + 1 ]);
3178               for ( int k = 0; k < 2; ++k )
3179               {
3180                 corner.SetCoord( 3, sP[ i*3 + 2 ]);
3181                 corner *= bi;
3182                 bb.Add( corner );
3183               }
3184             }
3185           }
3186           corner = isMax ? bb.CornerMax() : bb.CornerMin();
3187           pAxis.SetCoord( iDir+1, corner.Coord( iDir+1 ));
3188         }
3189         else
3190         {
3191           gp_XYZ pf = pFaces.XYZ() * bi;
3192           pAxis.SetCoord( iDir+1, pf.Coord( iDir+1 ) );
3193         }
3194       }
3195     } // loop on 3 axes
3196
3197     shapeBox.SetVoid();
3198     shapeBox.Add( pMin );
3199     shapeBox.Add( pMax );
3200
3201     return;
3202   }
3203
3204 } // namespace
3205
3206 //=============================================================================
3207 /*!
3208  * \brief Generates 3D structured Cartesian mesh in the internal part of
3209  * solid shapes and polyhedral volumes near the shape boundary.
3210  *  \param theMesh - mesh to fill in
3211  *  \param theShape - a compound of all SOLIDs to mesh
3212  *  \retval bool - true in case of success
3213  */
3214 //=============================================================================
3215
3216 bool StdMeshers_Cartesian_3D::Compute(SMESH_Mesh &         theMesh,
3217                                       const TopoDS_Shape & theShape)
3218 {
3219   // The algorithm generates the mesh in following steps:
3220
3221   // 1) Intersection of grid lines with the geometry boundary.
3222   // This step allows to find out if a given node of the initial grid is
3223   // inside or outside the geometry.
3224
3225   // 2) For each cell of the grid, check how many of it's nodes are outside
3226   // of the geometry boundary. Depending on a result of this check
3227   // - skip a cell, if all it's nodes are outside
3228   // - skip a cell, if it is too small according to the size threshold
3229   // - add a hexahedron in the mesh, if all nodes are inside
3230   // - add a polyhedron in the mesh, if some nodes are inside and some outside
3231
3232   _computeCanceled = false;
3233
3234   SMESH_MesherHelper helper( theMesh );
3235
3236   try
3237   {
3238     Grid grid;
3239     grid._helper = &helper;
3240
3241     vector< TopoDS_Shape > faceVec;
3242     {
3243       TopTools_MapOfShape faceMap;
3244       TopExp_Explorer fExp;
3245       for ( fExp.Init( theShape, TopAbs_FACE ); fExp.More(); fExp.Next() )
3246         if ( !faceMap.Add( fExp.Current() ))
3247           faceMap.Remove( fExp.Current() ); // remove a face shared by two solids
3248
3249       for ( fExp.ReInit(); fExp.More(); fExp.Next() )
3250         if ( faceMap.Contains( fExp.Current() ))
3251           faceVec.push_back( fExp.Current() );
3252     }
3253     vector<FaceGridIntersector> facesItersectors( faceVec.size() );
3254     map< TGeomID, vector< TGeomID > > edge2faceIDsMap;
3255     TopExp_Explorer eExp;
3256     Bnd_Box shapeBox;
3257     for ( int i = 0; i < faceVec.size(); ++i )
3258     {
3259       facesItersectors[i]._face   = TopoDS::Face    ( faceVec[i] );
3260       facesItersectors[i]._faceID = grid._shapes.Add( faceVec[i] );
3261       facesItersectors[i]._grid   = &grid;
3262       shapeBox.Add( facesItersectors[i].GetFaceBndBox() );
3263
3264       if ( _hyp->GetToAddEdges() )
3265       {
3266         helper.SetSubShape( faceVec[i] );
3267         for ( eExp.Init( faceVec[i], TopAbs_EDGE ); eExp.More(); eExp.Next() )
3268         {
3269           const TopoDS_Edge& edge = TopoDS::Edge( eExp.Current() );
3270           if ( !SMESH_Algo::isDegenerated( edge ) &&
3271                !helper.IsRealSeam( edge ))
3272             edge2faceIDsMap[ grid._shapes.Add( edge )].push_back( facesItersectors[i]._faceID );
3273         }
3274       }
3275     }
3276
3277     getExactBndBox( faceVec, _hyp->GetAxisDirs(), shapeBox );
3278
3279     vector<double> xCoords, yCoords, zCoords;
3280     _hyp->GetCoordinates( xCoords, yCoords, zCoords, shapeBox );
3281
3282     grid.SetCoordinates( xCoords, yCoords, zCoords, _hyp->GetAxisDirs(), shapeBox );
3283
3284     if ( _computeCanceled ) return false;
3285
3286 #ifdef WITH_TBB
3287     { // copy partner faces and curves of not thread-safe types
3288       set< const Standard_Transient* > tshapes;
3289       BRepBuilderAPI_Copy copier;
3290       for ( size_t i = 0; i < facesItersectors.size(); ++i )
3291       {
3292         if ( !facesItersectors[i].IsThreadSafe(tshapes) )
3293         {
3294           copier.Perform( facesItersectors[i]._face );
3295           facesItersectors[i]._face = TopoDS::Face( copier );
3296         }
3297       }
3298     }
3299     // Intersection of grid lines with the geometry boundary.
3300     tbb::parallel_for ( tbb::blocked_range<size_t>( 0, facesItersectors.size() ),
3301                         ParallelIntersector( facesItersectors ),
3302                         tbb::simple_partitioner());
3303 #else
3304     for ( size_t i = 0; i < facesItersectors.size(); ++i )
3305       facesItersectors[i].Intersect();
3306 #endif
3307
3308     // put interesection points onto the GridLine's; this is done after intersection
3309     // to avoid contention of facesItersectors for writing into the same GridLine
3310     // in case of parallel work of facesItersectors
3311     for ( size_t i = 0; i < facesItersectors.size(); ++i )
3312       facesItersectors[i].StoreIntersections();
3313
3314     TopExp_Explorer solidExp (theShape, TopAbs_SOLID);
3315     helper.SetSubShape( solidExp.Current() );
3316     helper.SetElementsOnShape( true );
3317
3318     if ( _computeCanceled ) return false;
3319
3320     // create nodes on the geometry
3321     grid.ComputeNodes(helper);
3322
3323     if ( _computeCanceled ) return false;
3324
3325     // create volume elements
3326     Hexahedron hex( _hyp->GetSizeThreshold(), &grid );
3327     int nbAdded = hex.MakeElements( helper, edge2faceIDsMap );
3328
3329     SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
3330     if ( nbAdded > 0 )
3331     {
3332       // make all SOLIDs computed
3333       if ( SMESHDS_SubMesh* sm1 = meshDS->MeshElements( solidExp.Current()) )
3334       {
3335         SMDS_ElemIteratorPtr volIt = sm1->GetElements();
3336         for ( ; solidExp.More() && volIt->more(); solidExp.Next() )
3337         {
3338           const SMDS_MeshElement* vol = volIt->next();
3339           sm1->RemoveElement( vol, /*isElemDeleted=*/false );
3340           meshDS->SetMeshElementOnShape( vol, solidExp.Current() );
3341         }
3342       }
3343       // make other sub-shapes computed
3344       setSubmeshesComputed( theMesh, theShape );
3345     }
3346
3347     // remove free nodes
3348     if ( SMESHDS_SubMesh * smDS = meshDS->MeshElements( helper.GetSubShapeID() ))
3349     {
3350       TIDSortedNodeSet nodesToRemove;
3351       // get intersection nodes
3352       for ( int iDir = 0; iDir < 3; ++iDir )
3353       {
3354         vector< GridLine >& lines = grid._lines[ iDir ];
3355         for ( size_t i = 0; i < lines.size(); ++i )
3356         {
3357           multiset< F_IntersectPoint >::iterator ip = lines[i]._intPoints.begin();
3358           for ( ; ip != lines[i]._intPoints.end(); ++ip )
3359             if ( ip->_node && ip->_node->NbInverseElements() == 0 )
3360               nodesToRemove.insert( nodesToRemove.end(), ip->_node );
3361         }
3362       }
3363       // get grid nodes
3364       for ( size_t i = 0; i < grid._nodes.size(); ++i )
3365         if ( grid._nodes[i] && grid._nodes[i]->NbInverseElements() == 0 )
3366           nodesToRemove.insert( nodesToRemove.end(), grid._nodes[i] );
3367
3368       // do remove
3369       TIDSortedNodeSet::iterator n = nodesToRemove.begin();
3370       for ( ; n != nodesToRemove.end(); ++n )
3371         meshDS->RemoveFreeNode( *n, smDS, /*fromGroups=*/false );
3372     }
3373
3374     return nbAdded;
3375
3376   }
3377   // SMESH_ComputeError is not caught at SMESH_submesh level for an unknown reason
3378   catch ( SMESH_ComputeError& e)
3379   {
3380     return error( SMESH_ComputeErrorPtr( new SMESH_ComputeError( e )));
3381   }
3382   return false;
3383 }
3384
3385 //=============================================================================
3386 /*!
3387  *  Evaluate
3388  */
3389 //=============================================================================
3390
3391 bool StdMeshers_Cartesian_3D::Evaluate(SMESH_Mesh &         theMesh,
3392                                        const TopoDS_Shape & theShape,
3393                                        MapShapeNbElems&     theResMap)
3394 {
3395   // TODO
3396 //   std::vector<int> aResVec(SMDSEntity_Last);
3397 //   for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
3398 //   if(IsQuadratic) {
3399 //     aResVec[SMDSEntity_Quad_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
3400 //     int nb1d_face0_int = ( nb2d_face0*4 - nb1d ) / 2;
3401 //     aResVec[SMDSEntity_Node] = nb0d_face0 * ( 2*nb2d/nb1d - 1 ) - nb1d_face0_int * nb2d/nb1d;
3402 //   }
3403 //   else {
3404 //     aResVec[SMDSEntity_Node] = nb0d_face0 * ( nb2d/nb1d - 1 );
3405 //     aResVec[SMDSEntity_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
3406 //   }
3407 //   SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
3408 //   aResMap.insert(std::make_pair(sm,aResVec));
3409
3410   return true;
3411 }
3412
3413 //=============================================================================
3414 namespace
3415 {
3416   /*!
3417    * \brief Event listener setting/unsetting _alwaysComputed flag to
3418    *        submeshes of inferior levels to prevent their computing
3419    */
3420   struct _EventListener : public SMESH_subMeshEventListener
3421   {
3422     string _algoName;
3423
3424     _EventListener(const string& algoName):
3425       SMESH_subMeshEventListener(/*isDeletable=*/true,"StdMeshers_Cartesian_3D::_EventListener"),
3426       _algoName(algoName)
3427     {}
3428     // --------------------------------------------------------------------------------
3429     // setting/unsetting _alwaysComputed flag to submeshes of inferior levels
3430     //
3431     static void setAlwaysComputed( const bool     isComputed,
3432                                    SMESH_subMesh* subMeshOfSolid)
3433     {
3434       SMESH_subMeshIteratorPtr smIt =
3435         subMeshOfSolid->getDependsOnIterator(/*includeSelf=*/false, /*complexShapeFirst=*/false);
3436       while ( smIt->more() )
3437       {
3438         SMESH_subMesh* sm = smIt->next();
3439         sm->SetIsAlwaysComputed( isComputed );
3440       }
3441       subMeshOfSolid->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
3442     }
3443
3444     // --------------------------------------------------------------------------------
3445     // unsetting _alwaysComputed flag if "Cartesian_3D" was removed
3446     //
3447     virtual void ProcessEvent(const int          event,
3448                               const int          eventType,
3449                               SMESH_subMesh*     subMeshOfSolid,
3450                               SMESH_subMeshEventListenerData* data,
3451                               const SMESH_Hypothesis*         hyp = 0)
3452     {
3453       if ( eventType == SMESH_subMesh::COMPUTE_EVENT )
3454       {
3455         setAlwaysComputed( subMeshOfSolid->GetComputeState() == SMESH_subMesh::COMPUTE_OK,
3456                            subMeshOfSolid );
3457       }
3458       else
3459       {
3460         SMESH_Algo* algo3D = subMeshOfSolid->GetAlgo();
3461         if ( !algo3D || _algoName != algo3D->GetName() )
3462           setAlwaysComputed( false, subMeshOfSolid );
3463       }
3464     }
3465
3466     // --------------------------------------------------------------------------------
3467     // set the event listener
3468     //
3469     static void SetOn( SMESH_subMesh* subMeshOfSolid, const string& algoName )
3470     {
3471       subMeshOfSolid->SetEventListener( new _EventListener( algoName ),
3472                                         /*data=*/0,
3473                                         subMeshOfSolid );
3474     }
3475
3476   }; // struct _EventListener
3477
3478 } // namespace
3479
3480 //================================================================================
3481 /*!
3482  * \brief Sets event listener to submeshes if necessary
3483  *  \param subMesh - submesh where algo is set
3484  * This method is called when a submesh gets HYP_OK algo_state.
3485  * After being set, event listener is notified on each event of a submesh.
3486  */
3487 //================================================================================
3488
3489 void StdMeshers_Cartesian_3D::SetEventListener(SMESH_subMesh* subMesh)
3490 {
3491   _EventListener::SetOn( subMesh, GetName() );
3492 }
3493
3494 //================================================================================
3495 /*!
3496  * \brief Set _alwaysComputed flag to submeshes of inferior levels to avoid their computing
3497  */
3498 //================================================================================
3499
3500 void StdMeshers_Cartesian_3D::setSubmeshesComputed(SMESH_Mesh&         theMesh,
3501                                                    const TopoDS_Shape& theShape)
3502 {
3503   for ( TopExp_Explorer soExp( theShape, TopAbs_SOLID ); soExp.More(); soExp.Next() )
3504     _EventListener::setAlwaysComputed( true, theMesh.GetSubMesh( soExp.Current() ));
3505 }
3506