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