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