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