Salome HOME
Merge from V6_main (04/10/2012)
[modules/smesh.git] / src / StdMeshers / StdMeshers_Cartesian_3D.cxx
1 // Copyright (C) 2007-2012  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
5 //
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License.
10 //
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
14 // Lesser General Public License for more details.
15 //
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
19 //
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22 //  File   : StdMeshers_Cartesian_3D.cxx
23 //  Module : SMESH
24 //
25 #include "StdMeshers_Cartesian_3D.hxx"
26
27 #include "SMDS_MeshNode.hxx"
28 #include "SMESH_Block.hxx"
29 #include "SMESH_Comment.hxx"
30 #include "SMESH_Mesh.hxx"
31 #include "SMESH_MesherHelper.hxx"
32 #include "SMESH_subMesh.hxx"
33 #include "SMESH_subMeshEventListener.hxx"
34 #include "StdMeshers_CartesianParameters3D.hxx"
35
36 #include "utilities.h"
37 #include "Utils_ExceptHandlers.hxx"
38 #include <Basics_OCCTVersion.hxx>
39
40 #include <BRepAdaptor_Surface.hxx>
41 #include <BRepBndLib.hxx>
42 #include <BRepBuilderAPI_Copy.hxx>
43 #include <BRepTools.hxx>
44 #include <BRep_Tool.hxx>
45 #include <Bnd_Box.hxx>
46 #include <ElSLib.hxx>
47 #include <Geom2d_BSplineCurve.hxx>
48 #include <Geom2d_BezierCurve.hxx>
49 #include <Geom2d_TrimmedCurve.hxx>
50 #include <Geom_BSplineCurve.hxx>
51 #include <Geom_BSplineSurface.hxx>
52 #include <Geom_BezierCurve.hxx>
53 #include <Geom_BezierSurface.hxx>
54 #include <Geom_RectangularTrimmedSurface.hxx>
55 #include <Geom_TrimmedCurve.hxx>
56 #include <IntAna_IntConicQuad.hxx>
57 #include <IntAna_IntLinTorus.hxx>
58 #include <IntAna_Quadric.hxx>
59 #include <IntCurveSurface_TransitionOnCurve.hxx>
60 #include <IntCurvesFace_Intersector.hxx>
61 #include <Poly_Triangulation.hxx>
62 #include <Precision.hxx>
63 #include <TopExp.hxx>
64 #include <TopExp_Explorer.hxx>
65 #include <TopLoc_Location.hxx>
66 #include <TopTools_MapIteratorOfMapOfShape.hxx>
67 #include <TopTools_MapOfShape.hxx>
68 #include <TopoDS.hxx>
69 #include <TopoDS_Face.hxx>
70 #include <TopoDS_TShape.hxx>
71 #include <gp_Cone.hxx>
72 #include <gp_Cylinder.hxx>
73 #include <gp_Lin.hxx>
74 #include <gp_Pln.hxx>
75 #include <gp_Pnt2d.hxx>
76 #include <gp_Sphere.hxx>
77 #include <gp_Torus.hxx>
78
79 //#undef WITH_TBB
80 #ifdef WITH_TBB
81 #include <tbb/parallel_for.h>
82 //#include <tbb/enumerable_thread_specific.h>
83 #endif
84
85 using namespace std;
86
87 //#define _MY_DEBUG_
88
89 #if OCC_VERSION_LARGE <= 0x06050300
90 // workaround it required only for OCCT6.5.3 and older (see OCC22809)
91 #define ELLIPSOLID_WORKAROUND
92 #endif
93
94 #ifdef ELLIPSOLID_WORKAROUND
95 #include <BRepIntCurveSurface_Inter.hxx>
96 #include <BRepTopAdaptor_TopolTool.hxx>
97 #include <BRepAdaptor_HSurface.hxx>
98 #endif
99
100 //=============================================================================
101 /*!
102  * Constructor
103  */
104 //=============================================================================
105
106 StdMeshers_Cartesian_3D::StdMeshers_Cartesian_3D(int hypId, int studyId, SMESH_Gen * gen)
107   :SMESH_3D_Algo(hypId, studyId, gen)
108 {
109   _name = "Cartesian_3D";
110   _shapeType = (1 << TopAbs_SOLID);       // 1 bit /shape type
111   _compatibleHypothesis.push_back("CartesianParameters3D");
112
113   _onlyUnaryInput = false;          // to mesh all SOLIDs at once
114   _requireDiscreteBoundary = false; // 2D mesh not needed
115   _supportSubmeshes = false;        // do not use any existing mesh
116 }
117
118 //=============================================================================
119 /*!
120  * Check presence of a hypothesis
121  */
122 //=============================================================================
123
124 bool StdMeshers_Cartesian_3D::CheckHypothesis (SMESH_Mesh&          aMesh,
125                                                const TopoDS_Shape&  aShape,
126                                                Hypothesis_Status&   aStatus)
127 {
128   aStatus = SMESH_Hypothesis::HYP_MISSING;
129
130   const list<const SMESHDS_Hypothesis*>& hyps = GetUsedHypothesis(aMesh, aShape);
131   list <const SMESHDS_Hypothesis* >::const_iterator h = hyps.begin();
132   if ( h == hyps.end())
133   {
134     return false;
135   }
136
137   for ( ; h != hyps.end(); ++h )
138   {
139     if (( _hyp = dynamic_cast<const StdMeshers_CartesianParameters3D*>( *h )))
140     {
141       aStatus = _hyp->IsDefined() ? HYP_OK : HYP_BAD_PARAMETER;
142       break;
143     }
144   }
145
146   return aStatus == HYP_OK;
147 }
148
149 namespace
150 {
151   //=============================================================================
152   // Definitions of internal utils
153   // --------------------------------------------------------------------------
154   enum Transition {
155     Trans_TANGENT = IntCurveSurface_Tangent,
156     Trans_IN      = IntCurveSurface_In,
157     Trans_OUT     = IntCurveSurface_Out,
158     Trans_APEX
159   };
160   // --------------------------------------------------------------------------
161   /*!
162    * \brief Data of intersection between a GridLine and a TopoDS_Face
163    */
164   struct IntersectionPoint
165   {
166     double                       _paramOnLine;
167     mutable Transition           _transition;
168     mutable const SMDS_MeshNode* _node;
169     mutable size_t               _indexOnLine;
170
171     IntersectionPoint(): _node(0) {}
172     bool operator< ( const IntersectionPoint& o ) const { return _paramOnLine < o._paramOnLine; }
173   };
174   // --------------------------------------------------------------------------
175   /*!
176    * \brief A line of the grid and its intersections with 2D geometry
177    */
178   struct GridLine
179   {
180     gp_Lin _line;
181     double _length; // line length
182     multiset< IntersectionPoint > _intPoints;
183
184     void RemoveExcessIntPoints( const double tol );
185     bool GetIsOutBefore( multiset< IntersectionPoint >::iterator ip, bool prevIsOut );
186   };
187   // --------------------------------------------------------------------------
188   /*!
189    * \brief Iterator on the parallel grid lines of one direction
190    */
191   struct LineIndexer
192   {
193     size_t _size  [3];
194     size_t _curInd[3];
195     size_t _iVar1, _iVar2, _iConst;
196     string _name1, _name2, _nameConst;
197     LineIndexer() {}
198     LineIndexer( size_t sz1, size_t sz2, size_t sz3,
199                  size_t iv1, size_t iv2, size_t iConst,
200                  const string& nv1, const string& nv2, const string& nConst )
201     {
202       _size[0] = sz1; _size[1] = sz2; _size[2] = sz3;
203       _curInd[0] = _curInd[1] = _curInd[2] = 0;
204       _iVar1 = iv1; _iVar2 = iv2; _iConst = iConst; 
205       _name1 = nv1; _name2 = nv2; _nameConst = nConst;
206     }
207
208     size_t I() const { return _curInd[0]; }
209     size_t J() const { return _curInd[1]; }
210     size_t K() const { return _curInd[2]; }
211     void SetIJK( size_t i, size_t j, size_t k )
212     {
213       _curInd[0] = i; _curInd[1] = j; _curInd[2] = k;
214     }
215     void operator++()
216     {
217       if ( ++_curInd[_iVar1] == _size[_iVar1] )
218         _curInd[_iVar1] = 0, ++_curInd[_iVar2];
219     }
220     bool More() const { return _curInd[_iVar2] < _size[_iVar2]; }
221     size_t LineIndex   () const { return _curInd[_iVar1] + _curInd[_iVar2]* _size[_iVar1]; }
222     size_t LineIndex10 () const { return (_curInd[_iVar1] + 1 ) + _curInd[_iVar2]* _size[_iVar1]; }
223     size_t LineIndex01 () const { return _curInd[_iVar1] + (_curInd[_iVar2] + 1 )* _size[_iVar1]; }
224     size_t LineIndex11 () const { return (_curInd[_iVar1] + 1 ) + (_curInd[_iVar2] + 1 )* _size[_iVar1]; }
225     void SetIndexOnLine (size_t i)  { _curInd[ _iConst ] = i; }
226     size_t NbLines() const { return _size[_iVar1] * _size[_iVar2]; }
227   };
228   // --------------------------------------------------------------------------
229   /*!
230    * \brief Container of GridLine's
231    */
232   struct Grid
233   {
234     vector< double >   _coords[3]; // coordinates of grid nodes
235     vector< GridLine > _lines [3]; //  in 3 directions
236     double             _tol, _minCellSize;
237
238     vector< const SMDS_MeshNode* > _nodes; // mesh nodes at grid nodes
239     vector< bool >                 _isBndNode; // is mesh node at intersection with geometry
240
241     size_t CellIndex( size_t i, size_t j, size_t k ) const
242     {
243       return i + j*(_coords[0].size()-1) + k*(_coords[0].size()-1)*(_coords[1].size()-1);
244     }
245     size_t NodeIndex( size_t i, size_t j, size_t k ) const
246     {
247       return i + j*_coords[0].size() + k*_coords[0].size()*_coords[1].size();
248     }
249     size_t NodeIndexDX() const { return 1; }
250     size_t NodeIndexDY() const { return _coords[0].size(); }
251     size_t NodeIndexDZ() const { return _coords[0].size() * _coords[1].size(); }
252
253     LineIndexer GetLineIndexer(size_t iDir) const;
254
255     void SetCoordinates(const vector<double>& xCoords,
256                         const vector<double>& yCoords,
257                         const vector<double>& zCoords,
258                         const TopoDS_Shape&   shape );
259     void ComputeNodes(SMESH_MesherHelper& helper);
260   };
261 #ifdef ELLIPSOLID_WORKAROUND
262   // --------------------------------------------------------------------------
263   /*!
264    * \brief struct temporary replacing IntCurvesFace_Intersector until
265    *        OCCT bug 0022809 is fixed
266    *        http://tracker.dev.opencascade.org/view.php?id=22809
267    */
268   struct TMP_IntCurvesFace_Intersector
269   {
270     BRepAdaptor_Surface                       _surf;
271     double                                    _tol;
272     BRepIntCurveSurface_Inter                 _intcs;
273     vector<IntCurveSurface_IntersectionPoint> _points;
274     BRepTopAdaptor_TopolTool                  _clsf;
275
276     TMP_IntCurvesFace_Intersector(const TopoDS_Face& face, const double tol)
277       :_surf( face ), _tol( tol ), _clsf( new BRepAdaptor_HSurface(_surf) ) {}
278     Bnd_Box Bounding() const { Bnd_Box b; BRepBndLib::Add (_surf.Face(), b); return b; }
279     void Perform( const gp_Lin& line, const double w0, const double w1 )
280     {
281       _points.clear();
282       for ( _intcs.Init( _surf.Face(), line, _tol ); _intcs.More(); _intcs.Next() )
283         if ( w0 <= _intcs.W() && _intcs.W() <= w1 )
284           _points.push_back( _intcs.Point() );
285     }
286     bool IsDone() const { return true; }
287     int  NbPnt()  const { return _points.size(); }
288     IntCurveSurface_TransitionOnCurve Transition( const int i ) const { return _points[ i-1 ].Transition(); }
289     double       WParameter( const int i ) const { return _points[ i-1 ].W(); }
290     TopAbs_State ClassifyUVPoint(const gp_Pnt2d& p) { return _clsf.Classify( p, _tol ); }
291   };
292 #define __IntCurvesFace_Intersector TMP_IntCurvesFace_Intersector
293 #else
294 #define __IntCurvesFace_Intersector IntCurvesFace_Intersector
295 #endif
296   // --------------------------------------------------------------------------
297   /*!
298    * \brief Intersector of TopoDS_Face with all GridLine's
299    */
300   struct FaceGridIntersector
301   {
302     TopoDS_Face _face;
303     Grid*       _grid;
304     Bnd_Box     _bndBox;
305     __IntCurvesFace_Intersector* _surfaceInt;
306     vector< std::pair< GridLine*, IntersectionPoint > > _intersections;
307
308     FaceGridIntersector(): _grid(0), _surfaceInt(0) {}
309     void Intersect();
310     bool IsInGrid(const Bnd_Box& gridBox);
311
312     void StoreIntersections()
313     {
314       for ( size_t i = 0; i < _intersections.size(); ++i )
315         _intersections[i].first->_intPoints.insert( _intersections[i].second );
316     }
317     const Bnd_Box& GetFaceBndBox()
318     {
319       GetCurveFaceIntersector();
320       return _bndBox;
321     }
322     __IntCurvesFace_Intersector* GetCurveFaceIntersector()
323     {
324       if ( !_surfaceInt )
325       {
326         _surfaceInt = new __IntCurvesFace_Intersector( _face, Precision::PConfusion() );
327         _bndBox     = _surfaceInt->Bounding();
328         if ( _bndBox.IsVoid() )
329           BRepBndLib::Add (_face, _bndBox);
330       }
331       return _surfaceInt;
332     }
333     bool IsThreadSafe(set< const Standard_Transient* >& noSafeTShapes) const;
334   };
335   // --------------------------------------------------------------------------
336   /*!
337    * \brief Intersector of a surface with a GridLine
338    */
339   struct FaceLineIntersector
340   {
341     double      _tol;
342     double      _u, _v, _w; // params on the face and the line
343     Transition  _transition; // transition of at intersection (see IntCurveSurface.cdl)
344     Transition  _transIn, _transOut; // IN and OUT transitions depending of face orientation
345
346     gp_Pln      _plane;
347     gp_Cylinder _cylinder;
348     gp_Cone     _cone;
349     gp_Sphere   _sphere;
350     gp_Torus    _torus;
351     __IntCurvesFace_Intersector* _surfaceInt;
352
353     vector< IntersectionPoint > _intPoints;
354
355     void IntersectWithPlane   (const GridLine& gridLine);
356     void IntersectWithCylinder(const GridLine& gridLine);
357     void IntersectWithCone    (const GridLine& gridLine);
358     void IntersectWithSphere  (const GridLine& gridLine);
359     void IntersectWithTorus   (const GridLine& gridLine);
360     void IntersectWithSurface (const GridLine& gridLine);
361
362     bool UVIsOnFace() const;
363     void addIntPoint(const bool toClassify=true);
364     bool isParamOnLineOK( const double linLength )
365     {
366       return -_tol < _w && _w < linLength + _tol;
367     }
368     FaceLineIntersector():_surfaceInt(0) {}
369     ~FaceLineIntersector() { if (_surfaceInt ) delete _surfaceInt; _surfaceInt = 0; }
370   };
371   // --------------------------------------------------------------------------
372   /*!
373    * \brief Class representing topology of the hexahedron and creating a mesh
374    *        volume basing on analysis of hexahedron intersection with geometry
375    */
376   class Hexahedron
377   {
378     // --------------------------------------------------------------------------------
379     struct _Face;
380     struct _Link;
381     // --------------------------------------------------------------------------------
382     struct _Node //!< node either at a hexahedron corner or at GridLine intersection
383     {
384       const SMDS_MeshNode*     _node; // mesh node at hexahedron corner
385       const IntersectionPoint* _intPoint;
386
387       _Node(const SMDS_MeshNode* n=0, const IntersectionPoint* ip=0):_node(n), _intPoint(ip) {} 
388       const SMDS_MeshNode* Node() const { return _intPoint ? _intPoint->_node : _node; }
389       //bool IsCorner() const { return _node; }
390     };
391     // --------------------------------------------------------------------------------
392     struct _Link // link connecting two _Node's
393     {
394       _Node* _nodes[2];
395       vector< _Node>  _intNodes; // _Node's at GridLine intersections
396       vector< _Link > _splits;
397       vector< _Face*> _faces;
398     };
399     // --------------------------------------------------------------------------------
400     struct _OrientedLink
401     {
402       _Link* _link;
403       bool   _reverse;
404       _OrientedLink( _Link* link=0, bool reverse=false ): _link(link), _reverse(reverse) {}
405       void Reverse() { _reverse = !_reverse; }
406       int NbResultLinks() const { return _link->_splits.size(); }
407       _OrientedLink ResultLink(int i) const
408       {
409         return _OrientedLink(&_link->_splits[_reverse ? NbResultLinks()-i-1 : i],_reverse);
410       }
411       _Node* FirstNode() const { return _link->_nodes[ _reverse ]; }
412       _Node* LastNode() const { return _link->_nodes[ !_reverse ]; }
413     };
414     // --------------------------------------------------------------------------------
415     struct _Face
416     {
417       vector< _OrientedLink > _links;
418       vector< _Link >         _polyLinks; // links added to close a polygonal face
419     };
420     // --------------------------------------------------------------------------------
421     struct _volumeDef // holder of nodes of a volume mesh element
422     {
423       vector< const SMDS_MeshNode* > _nodes;
424       vector< int >                  _quantities;
425       typedef boost::shared_ptr<_volumeDef> Ptr;
426       void set( const vector< const SMDS_MeshNode* >& nodes,
427                 const vector< int > quant = vector< int >() )
428       { _nodes = nodes; _quantities = quant; }
429       // static Ptr New( const vector< const SMDS_MeshNode* >& nodes,
430       //                 const vector< int > quant = vector< int >() )
431       // {
432       //   _volumeDef* def = new _volumeDef;
433       //   def->_nodes = nodes;
434       //   def->_quantities = quant;
435       //   return Ptr( def );
436       // }
437     };
438
439     // topology of a hexahedron
440     int   _nodeShift[8];
441     _Node _hexNodes[8];
442     _Link _hexLinks[12];
443     _Face _hexQuads[6];
444
445     // faces resulted from hexahedron intersection
446     vector< _Face > _polygons;
447
448     // computed volume elements
449     //vector< _volumeDef::Ptr > _volumeDefs;
450     _volumeDef _volumeDefs;
451
452     Grid*       _grid;
453     double      _sizeThreshold, _sideLength[3];
454     int         _nbCornerNodes, _nbIntNodes, _nbBndNodes;
455     int         _origNodeInd; // index of _hexNodes[0] node within the _grid
456     size_t      _i,_j,_k;
457
458   public:
459     Hexahedron(const double sizeThreshold, Grid* grid);
460     int MakeElements(SMESH_MesherHelper& helper);
461     void ComputeElements();
462     void Init() { init( _i, _j, _k ); }
463
464   private:
465     Hexahedron(const Hexahedron& other );
466     void init( size_t i, size_t j, size_t k );
467     void init( size_t i );
468     int  addElements(SMESH_MesherHelper& helper);
469     bool isInHole() const;
470     bool checkPolyhedronSize() const;
471     bool addHexa ();
472     bool addTetra();
473     bool addPenta();
474     bool addPyra ();
475   };
476  
477 #ifdef WITH_TBB
478   // --------------------------------------------------------------------------
479   /*!
480    * \brief Hexahedron computing volumes in one thread
481    */
482   struct ParallelHexahedron
483   {
484     vector< Hexahedron* >& _hexVec;
485     vector<int>&           _index;
486     ParallelHexahedron( vector< Hexahedron* >& hv, vector<int>& ind): _hexVec(hv), _index(ind) {}
487     void operator() ( const tbb::blocked_range<size_t>& r ) const
488     {
489       for ( size_t i = r.begin(); i != r.end(); ++i )
490         if ( Hexahedron* hex = _hexVec[ _index[i]] )
491           hex->ComputeElements();
492     }
493   };
494   // --------------------------------------------------------------------------
495   /*!
496    * \brief Structure intersecting certain nb of faces with GridLine's in one thread
497    */
498   struct ParallelIntersector
499   {
500     vector< FaceGridIntersector >& _faceVec;
501     ParallelIntersector( vector< FaceGridIntersector >& faceVec): _faceVec(faceVec){}
502     void operator() ( const tbb::blocked_range<size_t>& r ) const
503     {
504       for ( size_t i = r.begin(); i != r.end(); ++i )
505         _faceVec[i].Intersect();
506     }
507   };
508
509 #endif
510   //=============================================================================
511   // Implementation of internal utils
512   //=============================================================================
513   /*
514    * Remove coincident intersection points
515    */
516   void GridLine::RemoveExcessIntPoints( const double tol )
517   {
518     if ( _intPoints.size() < 2 ) return;
519
520     set< Transition > tranSet;
521     multiset< IntersectionPoint >::iterator ip1, ip2 = _intPoints.begin();
522     while ( ip2 != _intPoints.end() )
523     {
524       tranSet.clear();
525       ip1 = ip2++;
526       while ( ip2->_paramOnLine - ip1->_paramOnLine <= tol  && ip2 != _intPoints.end())
527       {
528         tranSet.insert( ip1->_transition );
529         tranSet.insert( ip2->_transition );
530         _intPoints.erase( ip1 );
531         ip1 = ip2++;
532       }
533       if ( tranSet.size() > 1 ) // points with different transition coincide
534       {
535         bool isIN  = tranSet.count( Trans_IN );
536         bool isOUT = tranSet.count( Trans_OUT );
537         if ( isIN && isOUT )
538           (*ip1)._transition = Trans_TANGENT;
539         else
540           (*ip1)._transition = isIN ? Trans_IN : Trans_OUT;
541       }
542     }
543   }
544   //================================================================================
545   /*
546    * Return "is OUT" state for nodes before the given intersection point
547    */
548   bool GridLine::GetIsOutBefore( multiset< IntersectionPoint >::iterator ip, bool prevIsOut )
549   {
550     if ( ip->_transition == Trans_IN )
551       return true;
552     if ( ip->_transition == Trans_OUT )
553       return false;
554     if ( ip->_transition == Trans_APEX )
555     {
556       // singularity point (apex of a cone)
557       if ( _intPoints.size() == 1 || ip == _intPoints.begin() )
558         return true;
559       multiset< IntersectionPoint >::iterator ipBef = ip, ipAft = ++ip;
560       if ( ipAft == _intPoints.end() )
561         return false;
562       --ipBef;
563       if ( ipBef->_transition != ipAft->_transition )
564         return ( ipBef->_transition == Trans_OUT );
565       return ( ipBef->_transition != Trans_OUT );
566     }
567     return prevIsOut; // _transition == Trans_TANGENT
568   }
569   //================================================================================
570   /*
571    * Return an iterator on GridLine's in a given direction
572    */
573   LineIndexer Grid::GetLineIndexer(size_t iDir) const
574   {
575     const size_t indices[] = { 1,2,0, 0,2,1, 0,1,2 };
576     const string s[] = { "X", "Y", "Z" };
577     LineIndexer li( _coords[0].size(),  _coords[1].size(),    _coords[2].size(),
578                     indices[iDir*3],    indices[iDir*3+1],    indices[iDir*3+2],
579                     s[indices[iDir*3]], s[indices[iDir*3+1]], s[indices[iDir*3+2]]);
580     return li;
581   }
582   //=============================================================================
583   /*
584    * Creates GridLine's of the grid
585    */
586   void Grid::SetCoordinates(const vector<double>& xCoords,
587                             const vector<double>& yCoords,
588                             const vector<double>& zCoords,
589                             const TopoDS_Shape&   shape)
590   {
591     _coords[0] = xCoords;
592     _coords[1] = yCoords;
593     _coords[2] = zCoords;
594
595     // compute tolerance
596     _minCellSize = Precision::Infinite();
597     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
598     {
599       for ( size_t i = 1; i < _coords[ iDir ].size(); ++i )
600       {
601         double cellLen = _coords[ iDir ][ i ] - _coords[ iDir ][ i-1 ];
602         if ( cellLen < _minCellSize )
603           _minCellSize = cellLen;
604       }
605     }
606     if ( _minCellSize < Precision::Confusion() )
607       throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
608                                 SMESH_Comment("Too small cell size: ") << _tol );
609     _tol = _minCellSize / 1000.;
610
611     // attune grid extremities to shape bounding box computed by vertices
612     Bnd_Box shapeBox;
613     for ( TopExp_Explorer vExp( shape, TopAbs_VERTEX ); vExp.More(); vExp.Next() )
614       shapeBox.Add( BRep_Tool::Pnt( TopoDS::Vertex( vExp.Current() )));
615     
616     double sP[6]; // aXmin, aYmin, aZmin, aXmax, aYmax, aZmax
617     shapeBox.Get(sP[0],sP[1],sP[2],sP[3],sP[4],sP[5]);
618     double* cP[6] = { &_coords[0].front(), &_coords[1].front(), &_coords[2].front(),
619                       &_coords[0].back(),  &_coords[1].back(),  &_coords[2].back() };
620     for ( int i = 0; i < 6; ++i )
621       if ( fabs( sP[i] - *cP[i] ) < _tol )
622         *cP[i] = sP[i] + _tol/1000. * ( i < 3 ? +1 : -1 );
623
624     // create lines
625     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
626     {
627       LineIndexer li = GetLineIndexer( iDir );
628       _lines[iDir].resize( li.NbLines() );
629       double len = _coords[ iDir ].back() - _coords[iDir].front();
630       gp_Vec dir( iDir==0, iDir==1, iDir==2 );
631       for ( ; li.More(); ++li )
632       {
633         GridLine& gl = _lines[iDir][ li.LineIndex() ];
634         gl._line.SetLocation(gp_Pnt(_coords[0][li.I()], _coords[1][li.J()], _coords[2][li.K()])); 
635         gl._line.SetDirection( dir );
636         gl._length = len;
637       }
638     }
639   }
640   //================================================================================
641   /*
642    * Creates all nodes
643    */
644   void Grid::ComputeNodes(SMESH_MesherHelper& helper)
645   {
646     // state of each node of the grid relative to the geomerty
647     const size_t nbGridNodes = _coords[0].size() * _coords[1].size() * _coords[2].size();
648     vector< bool > isNodeOut( nbGridNodes, false );
649     _nodes.resize( nbGridNodes, 0 );
650     _isBndNode.resize( nbGridNodes, false );
651
652     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
653     {
654       LineIndexer li = GetLineIndexer( iDir );
655
656       // find out a shift of node index while walking along a GridLine in this direction
657       li.SetIndexOnLine( 0 );
658       size_t nIndex0 = NodeIndex( li.I(), li.J(), li.K() );
659       li.SetIndexOnLine( 1 );
660       const size_t nShift = NodeIndex( li.I(), li.J(), li.K() ) - nIndex0;
661       
662       const vector<double> & coords = _coords[ iDir ];
663       for ( ; li.More(); ++li ) // loop on lines in iDir
664       {
665         li.SetIndexOnLine( 0 );
666         nIndex0 = NodeIndex( li.I(), li.J(), li.K() );
667
668         GridLine& line = _lines[ iDir ][ li.LineIndex() ];
669         line.RemoveExcessIntPoints( _tol );
670         multiset< IntersectionPoint >& intPnts = _lines[ iDir ][ li.LineIndex() ]._intPoints;
671         multiset< IntersectionPoint >::iterator ip = intPnts.begin();
672
673         bool isOut = true;
674         const double* nodeCoord = & coords[0], *coord0 = nodeCoord, *coordEnd = coord0 + coords.size();
675         double nodeParam = 0;
676         for ( ; ip != intPnts.end(); ++ip )
677         {
678           // set OUT state or just skip IN nodes before ip
679           if ( nodeParam < ip->_paramOnLine - _tol )
680           {
681             isOut = line.GetIsOutBefore( ip, isOut );
682
683             while ( nodeParam < ip->_paramOnLine - _tol )
684             {
685               if ( isOut )
686                 isNodeOut[ nIndex0 + nShift * ( nodeCoord-coord0 ) ] = isOut;
687               if ( ++nodeCoord <  coordEnd )
688                 nodeParam = *nodeCoord - *coord0;
689               else
690                 break;
691             }
692             if ( nodeCoord == coordEnd ) break;
693           }
694           // create a mesh node on a GridLine at ip if it does not coincide with a grid node
695           if ( nodeParam > ip->_paramOnLine + _tol )
696           {
697             li.SetIndexOnLine( 0 );
698             double xyz[3] = { _coords[0][ li.I() ], _coords[1][ li.J() ], _coords[2][ li.K() ]};
699             xyz[ li._iConst ] += ip->_paramOnLine;
700             ip->_node = helper.AddNode( xyz[0], xyz[1], xyz[2] );
701             ip->_indexOnLine = nodeCoord-coord0-1;
702           }
703           // create a mesh node at ip concident with a grid node
704           else
705           {
706             int nodeIndex = nIndex0 + nShift * ( nodeCoord-coord0 );
707             if ( ! _nodes[ nodeIndex ] )
708             {
709               li.SetIndexOnLine( nodeCoord-coord0 );
710               double xyz[3] = { _coords[0][ li.I() ], _coords[1][ li.J() ], _coords[2][ li.K() ]};
711               _nodes[ nodeIndex ] = helper.AddNode( xyz[0], xyz[1], xyz[2] );
712               _isBndNode[ nodeIndex ] = true;
713             }
714             //ip->_node = _nodes[ nodeIndex ];
715             ip->_indexOnLine = nodeCoord-coord0;
716             if ( ++nodeCoord < coordEnd )
717               nodeParam = *nodeCoord - *coord0;
718           }
719         }
720         // set OUT state to nodes after the last ip
721         for ( ; nodeCoord < coordEnd; ++nodeCoord )
722           isNodeOut[ nIndex0 + nShift * ( nodeCoord-coord0 ) ] = true;
723       }
724     }
725
726     // Create mesh nodes at !OUT nodes of the grid
727
728     for ( size_t z = 0; z < _coords[2].size(); ++z )
729       for ( size_t y = 0; y < _coords[1].size(); ++y )
730         for ( size_t x = 0; x < _coords[0].size(); ++x )
731         {
732           size_t nodeIndex = NodeIndex( x, y, z );
733           if ( !isNodeOut[ nodeIndex ] && !_nodes[ nodeIndex] )
734             _nodes[ nodeIndex ] = helper.AddNode( _coords[0][x], _coords[1][y], _coords[2][z] );
735         }
736
737 #ifdef _MY_DEBUG_
738     // check validity of transitions
739     const char* trName[] = { "TANGENT", "IN", "OUT", "APEX" };
740     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
741     {
742       LineIndexer li = GetLineIndexer( iDir );
743       for ( ; li.More(); ++li )
744       {
745         multiset< IntersectionPoint >& intPnts = _lines[ iDir ][ li.LineIndex() ]._intPoints;
746         if ( intPnts.empty() ) continue;
747         if ( intPnts.size() == 1 )
748         {
749           if ( intPnts.begin()->_transition != Trans_TANGENT &&
750                intPnts.begin()->_transition != Trans_APEX )
751           throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
752                                     SMESH_Comment("Wrong SOLE transition of GridLine (")
753                                     << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
754                                     << ") along " << li._nameConst
755                                     << ": " << trName[ intPnts.begin()->_transition] );
756         }
757         else
758         {
759           if ( intPnts.begin()->_transition == Trans_OUT )
760             throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
761                                       SMESH_Comment("Wrong START transition of GridLine (")
762                                       << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
763                                       << ") along " << li._nameConst
764                                       << ": " << trName[ intPnts.begin()->_transition ]);
765           if ( intPnts.rbegin()->_transition == Trans_IN )
766             throw SMESH_ComputeError (COMPERR_ALGO_FAILED,
767                                       SMESH_Comment("Wrong END transition of GridLine (")
768                                       << li._curInd[li._iVar1] << ", " << li._curInd[li._iVar2]
769                                       << ") along " << li._nameConst
770                                     << ": " << trName[ intPnts.rbegin()->_transition ]);
771         }
772       }
773     }
774 #endif
775   }
776
777   //=============================================================================
778   /*
779    * Checks if the face is encosed by the grid
780    */
781   bool FaceGridIntersector::IsInGrid(const Bnd_Box& gridBox)
782   {
783     double x0,y0,z0, x1,y1,z1;
784     const Bnd_Box& faceBox = GetFaceBndBox();
785     faceBox.Get(x0,y0,z0, x1,y1,z1);
786
787     if ( !gridBox.IsOut( gp_Pnt( x0,y0,z0 )) &&
788          !gridBox.IsOut( gp_Pnt( x1,y1,z1 )))
789       return true;
790
791     double X0,Y0,Z0, X1,Y1,Z1;
792     gridBox.Get(X0,Y0,Z0, X1,Y1,Z1);
793     double faceP[6] = { x0,y0,z0, x1,y1,z1 };
794     double gridP[6] = { X0,Y0,Z0, X1,Y1,Z1 };
795     gp_Dir axes[3]  = { gp::DX(), gp::DY(), gp::DZ() };
796     for ( int iDir = 0; iDir < 6; ++iDir )
797     {
798       if ( iDir < 3  && gridP[ iDir ] <= faceP[ iDir ] ) continue;
799       if ( iDir >= 3 && gridP[ iDir ] >= faceP[ iDir ] ) continue;
800
801       // check if the face intersects a side of a gridBox
802
803       gp_Pnt p = iDir < 3 ? gp_Pnt( X0,Y0,Z0 ) : gp_Pnt( X1,Y1,Z1 );
804       gp_Ax1 norm( p, axes[ iDir % 3 ] );
805       if ( iDir < 3 ) norm.Reverse();
806
807       gp_XYZ O = norm.Location().XYZ(), N = norm.Direction().XYZ();
808
809       TopLoc_Location loc = _face.Location();
810       Handle(Poly_Triangulation) aPoly = BRep_Tool::Triangulation(_face,loc);
811       if ( !aPoly.IsNull() )
812       {
813         if ( !loc.IsIdentity() )
814         {
815           norm.Transform( loc.Transformation().Inverted() );
816           O = norm.Location().XYZ(), N = norm.Direction().XYZ();
817         }
818         const double deflection = aPoly->Deflection();
819
820         const TColgp_Array1OfPnt& nodes = aPoly->Nodes();
821         for ( int i = nodes.Lower(); i <= nodes.Upper(); ++i )
822           if (( nodes( i ).XYZ() - O ) * N > _grid->_tol + deflection )
823             return false;
824       }
825       else
826       {
827         BRepAdaptor_Surface surf( _face );
828         double u0, u1, v0, v1, du, dv, u, v;
829         BRepTools::UVBounds( _face, u0, u1, v0, v1);
830         if ( surf.GetType() == GeomAbs_Plane ) {
831           du = u1 - u0, dv = v1 - v0;
832         }
833         else {
834           du = surf.UResolution( _grid->_minCellSize / 10. );
835           dv = surf.VResolution( _grid->_minCellSize / 10. );
836         }
837         for ( u = u0, v = v0; u <= u1 && v <= v1; u += du, v += dv )
838         {
839           gp_Pnt p = surf.Value( u, v );
840           if (( p.XYZ() - O ) * N > _grid->_tol )
841           {
842             TopAbs_State state = GetCurveFaceIntersector()->ClassifyUVPoint(gp_Pnt2d( u, v ));
843             if ( state == TopAbs_IN || state == TopAbs_ON )
844               return false;
845           }
846         }
847       }
848     }
849     return true;
850   }
851   //=============================================================================
852   /*
853    * Intersects TopoDS_Face with all GridLine's
854    */
855   void FaceGridIntersector::Intersect()
856   {
857     FaceLineIntersector intersector;
858     intersector._surfaceInt = GetCurveFaceIntersector();
859     intersector._tol        = _grid->_tol;
860     intersector._transOut   = _face.Orientation() == TopAbs_REVERSED ? Trans_IN : Trans_OUT;
861     intersector._transIn    = _face.Orientation() == TopAbs_REVERSED ? Trans_OUT : Trans_IN;
862
863     typedef void (FaceLineIntersector::* PIntFun )(const GridLine& gridLine);
864     PIntFun interFunction;
865
866     BRepAdaptor_Surface surf( _face );
867     switch ( surf.GetType() ) {
868     case GeomAbs_Plane:
869       intersector._plane = surf.Plane();
870       interFunction = &FaceLineIntersector::IntersectWithPlane;
871       break;
872     case GeomAbs_Cylinder:
873       intersector._cylinder = surf.Cylinder();
874       interFunction = &FaceLineIntersector::IntersectWithCylinder;
875       break;
876     case GeomAbs_Cone:
877       intersector._cone = surf.Cone();
878       interFunction = &FaceLineIntersector::IntersectWithCone;
879       break;
880     case GeomAbs_Sphere:
881       intersector._sphere = surf.Sphere();
882       interFunction = &FaceLineIntersector::IntersectWithSphere;
883       break;
884     case GeomAbs_Torus:
885       intersector._torus = surf.Torus();
886       interFunction = &FaceLineIntersector::IntersectWithTorus;
887       break;
888     default:
889       interFunction = &FaceLineIntersector::IntersectWithSurface;
890     }
891
892     _intersections.clear();
893     for ( int iDir = 0; iDir < 3; ++iDir ) // loop on 3 line directions
894     {
895       if ( surf.GetType() == GeomAbs_Plane )
896       {
897         // check if all lines in this direction are parallel to a plane
898         if ( intersector._plane.Axis().IsNormal( _grid->_lines[iDir][0]._line.Position(),
899                                                  Precision::Angular()))
900           continue;
901         // find out a transition, that is the same for all lines of a direction
902         gp_Dir plnNorm = intersector._plane.Axis().Direction();
903         gp_Dir lineDir = _grid->_lines[iDir][0]._line.Direction();
904         intersector._transition =
905           ( plnNorm * lineDir < 0 ) ? intersector._transIn : intersector._transOut;
906       }
907       if ( surf.GetType() == GeomAbs_Cylinder )
908       {
909         // check if all lines in this direction are parallel to a cylinder
910         if ( intersector._cylinder.Axis().IsParallel( _grid->_lines[iDir][0]._line.Position(),
911                                                       Precision::Angular()))
912           continue;
913       }
914
915       // intersect the grid lines with the face
916       for ( size_t iL = 0; iL < _grid->_lines[iDir].size(); ++iL )
917       {
918         GridLine& gridLine = _grid->_lines[iDir][iL];
919         if ( _bndBox.IsOut( gridLine._line )) continue;
920
921         intersector._intPoints.clear();
922         (intersector.*interFunction)( gridLine );
923         for ( size_t i = 0; i < intersector._intPoints.size(); ++i )
924           _intersections.push_back( make_pair( &gridLine, intersector._intPoints[i] ));
925       }
926     }
927   }
928   //================================================================================
929   /*
930    * Return true if (_u,_v) is on the face
931    */
932   bool FaceLineIntersector::UVIsOnFace() const
933   {
934     TopAbs_State state = _surfaceInt->ClassifyUVPoint(gp_Pnt2d( _u,_v ));
935     return ( state == TopAbs_IN || state == TopAbs_ON );
936   }
937   //================================================================================
938   /*
939    * Store an intersection if it is IN or ON the face
940    */
941   void FaceLineIntersector::addIntPoint(const bool toClassify)
942   {
943     if ( !toClassify || UVIsOnFace() )
944     {
945       IntersectionPoint p;
946       p._paramOnLine = _w;
947       p._transition  = _transition;
948       _intPoints.push_back( p );
949     }
950   }
951   //================================================================================
952   /*
953    * Intersect a line with a plane
954    */
955   void FaceLineIntersector::IntersectWithPlane   (const GridLine& gridLine)
956   {
957     IntAna_IntConicQuad linPlane( gridLine._line, _plane, Precision::Angular());
958     _w = linPlane.ParamOnConic(1);
959     if ( isParamOnLineOK( gridLine._length ))
960     {
961       ElSLib::Parameters(_plane, linPlane.Point(1) ,_u,_v);
962       addIntPoint();
963     }
964   }
965   //================================================================================
966   /*
967    * Intersect a line with a cylinder
968    */
969   void FaceLineIntersector::IntersectWithCylinder(const GridLine& gridLine)
970   {
971     IntAna_IntConicQuad linCylinder( gridLine._line,_cylinder);
972     if ( linCylinder.IsDone() && linCylinder.NbPoints() > 0 )
973     {
974       _w = linCylinder.ParamOnConic(1);
975       if ( linCylinder.NbPoints() == 1 )
976         _transition = Trans_TANGENT;
977       else
978         _transition = _w < linCylinder.ParamOnConic(2) ? _transIn : _transOut;
979       if ( isParamOnLineOK( gridLine._length ))
980       {
981         ElSLib::Parameters(_cylinder, linCylinder.Point(1) ,_u,_v);
982         addIntPoint();
983       }
984       if ( linCylinder.NbPoints() > 1 )
985       {
986         _w = linCylinder.ParamOnConic(2);
987         if ( isParamOnLineOK( gridLine._length ))
988         {
989           ElSLib::Parameters(_cylinder, linCylinder.Point(2) ,_u,_v);
990           _transition = ( _transition == Trans_OUT ) ? Trans_IN : Trans_OUT;
991           addIntPoint();
992         }
993       }
994     }
995   }
996   //================================================================================
997   /*
998    * Intersect a line with a cone
999    */
1000   void FaceLineIntersector::IntersectWithCone (const GridLine& gridLine)
1001   {
1002     IntAna_IntConicQuad linCone(gridLine._line,_cone);
1003     if ( !linCone.IsDone() ) return;
1004     gp_Pnt P;
1005     gp_Vec du, dv, norm;
1006     for ( int i = 1; i <= linCone.NbPoints(); ++i )
1007     {
1008       _w = linCone.ParamOnConic( i );
1009       if ( !isParamOnLineOK( gridLine._length )) continue;
1010       ElSLib::Parameters(_cone, linCone.Point(i) ,_u,_v);
1011       if ( UVIsOnFace() )
1012       {
1013         ElSLib::D1( _u, _v, _cone, P, du, dv );
1014         norm = du ^ dv;
1015         double normSize2 = norm.SquareMagnitude();
1016         if ( normSize2 > Precision::Angular() * Precision::Angular() )
1017         {
1018           double cos = norm.XYZ() * gridLine._line.Direction().XYZ();
1019           cos /= sqrt( normSize2 );
1020           if ( cos < -Precision::Angular() )
1021             _transition = _transIn;
1022           else if ( cos > Precision::Angular() )
1023             _transition = _transOut;
1024           else
1025             _transition = Trans_TANGENT;
1026         }
1027         else
1028         {
1029           _transition = Trans_APEX;
1030         }
1031         addIntPoint( /*toClassify=*/false);
1032       }
1033     }
1034   }
1035   //================================================================================
1036   /*
1037    * Intersect a line with a sphere
1038    */
1039   void FaceLineIntersector::IntersectWithSphere  (const GridLine& gridLine)
1040   {
1041     IntAna_IntConicQuad linSphere(gridLine._line,_sphere);
1042     if ( linSphere.IsDone() && linSphere.NbPoints() > 0 )
1043     {
1044       _w = linSphere.ParamOnConic(1);
1045       if ( linSphere.NbPoints() == 1 )
1046         _transition = Trans_TANGENT;
1047       else
1048         _transition = _w < linSphere.ParamOnConic(2) ? _transIn : _transOut;
1049       if ( isParamOnLineOK( gridLine._length ))
1050       {
1051         ElSLib::Parameters(_sphere, linSphere.Point(1) ,_u,_v);
1052         addIntPoint();
1053       }
1054       if ( linSphere.NbPoints() > 1 )
1055       {
1056         _w = linSphere.ParamOnConic(2);
1057         if ( isParamOnLineOK( gridLine._length ))
1058         {
1059           ElSLib::Parameters(_sphere, linSphere.Point(2) ,_u,_v);
1060           _transition = ( _transition == Trans_OUT ) ? Trans_IN : Trans_OUT;
1061           addIntPoint();
1062         }
1063       }
1064     }
1065   }
1066   //================================================================================
1067   /*
1068    * Intersect a line with a torus
1069    */
1070   void FaceLineIntersector::IntersectWithTorus   (const GridLine& gridLine)
1071   {
1072     IntAna_IntLinTorus linTorus(gridLine._line,_torus);
1073     if ( !linTorus.IsDone()) return;
1074     gp_Pnt P;
1075     gp_Vec du, dv, norm;
1076     for ( int i = 1; i <= linTorus.NbPoints(); ++i )
1077     {
1078       _w = linTorus.ParamOnLine( i );
1079       if ( !isParamOnLineOK( gridLine._length )) continue;
1080       linTorus.ParamOnTorus( i, _u,_v );
1081       if ( UVIsOnFace() )
1082       {
1083         ElSLib::D1( _u, _v, _torus, P, du, dv );
1084         norm = du ^ dv;
1085         double normSize = norm.Magnitude();
1086         double cos = norm.XYZ() * gridLine._line.Direction().XYZ();
1087         cos /= normSize;
1088         if ( cos < -Precision::Angular() )
1089           _transition = _transIn;
1090         else if ( cos > Precision::Angular() )
1091           _transition = _transOut;
1092         else
1093           _transition = Trans_TANGENT;
1094         addIntPoint( /*toClassify=*/false);
1095       }
1096     }
1097   }
1098   //================================================================================
1099   /*
1100    * Intersect a line with a non-analytical surface
1101    */
1102   void FaceLineIntersector::IntersectWithSurface (const GridLine& gridLine)
1103   {
1104     _surfaceInt->Perform( gridLine._line, 0.0, gridLine._length );
1105     if ( !_surfaceInt->IsDone() ) return;
1106     for ( int i = 1; i <= _surfaceInt->NbPnt(); ++i )
1107     {
1108       _transition = Transition( _surfaceInt->Transition( i ) );
1109       _w = _surfaceInt->WParameter( i );
1110       addIntPoint(/*toClassify=*/false);
1111     }
1112   }
1113   //================================================================================
1114   /*
1115    * check if its face can be safely intersected in a thread
1116    */
1117   bool FaceGridIntersector::IsThreadSafe(set< const Standard_Transient* >& noSafeTShapes) const
1118   {
1119     bool isSafe = true;
1120
1121     // check surface
1122     TopLoc_Location loc;
1123     Handle(Geom_Surface) surf = BRep_Tool::Surface( _face, loc );
1124     Handle(Geom_RectangularTrimmedSurface) ts =
1125       Handle(Geom_RectangularTrimmedSurface)::DownCast( surf );
1126     while( !ts.IsNull() ) {
1127       surf = ts->BasisSurface();
1128       ts = Handle(Geom_RectangularTrimmedSurface)::DownCast(surf);
1129     }
1130     if ( surf->IsKind( STANDARD_TYPE(Geom_BSplineSurface )) ||
1131          surf->IsKind( STANDARD_TYPE(Geom_BezierSurface )))
1132       if ( !noSafeTShapes.insert((const Standard_Transient*) _face.TShape() ).second )
1133         isSafe = false;
1134
1135     double f, l;
1136     TopExp_Explorer exp( _face, TopAbs_EDGE );
1137     for ( ; exp.More(); exp.Next() )
1138     {
1139       bool edgeIsSafe = true;
1140       const TopoDS_Edge& e = TopoDS::Edge( exp.Current() );
1141       // check 3d curve
1142       {
1143         Handle(Geom_Curve) c = BRep_Tool::Curve( e, loc, f, l);
1144         if ( !c.IsNull() )
1145         {
1146           Handle(Geom_TrimmedCurve) tc = Handle(Geom_TrimmedCurve)::DownCast(c);
1147           while( !tc.IsNull() ) {
1148             c = tc->BasisCurve();
1149             tc = Handle(Geom_TrimmedCurve)::DownCast(c);
1150           }
1151           if ( c->IsKind( STANDARD_TYPE(Geom_BSplineCurve )) ||
1152                c->IsKind( STANDARD_TYPE(Geom_BezierCurve )))
1153             edgeIsSafe = false;
1154         }
1155       }
1156       // check 2d curve
1157       if ( edgeIsSafe )
1158       {
1159         Handle(Geom2d_Curve) c2 = BRep_Tool::CurveOnSurface( e, surf, loc, f, l);
1160         if ( !c2.IsNull() )
1161         {
1162           Handle(Geom2d_TrimmedCurve) tc = Handle(Geom2d_TrimmedCurve)::DownCast(c2);
1163           while( !tc.IsNull() ) {
1164             c2 = tc->BasisCurve();
1165             tc = Handle(Geom2d_TrimmedCurve)::DownCast(c2);
1166           }
1167           if ( c2->IsKind( STANDARD_TYPE(Geom2d_BSplineCurve )) ||
1168                c2->IsKind( STANDARD_TYPE(Geom2d_BezierCurve )))
1169             edgeIsSafe = false;
1170         }
1171       }
1172       if ( !edgeIsSafe && !noSafeTShapes.insert((const Standard_Transient*) e.TShape() ).second )
1173         isSafe = false;
1174     }
1175     return isSafe;
1176   }
1177   //================================================================================
1178   /*!
1179    * \brief Creates topology of the hexahedron
1180    */
1181   Hexahedron::Hexahedron(const double sizeThreshold, Grid* grid)
1182     : _grid( grid ), _sizeThreshold( sizeThreshold ), _nbIntNodes(0)
1183   {
1184     _polygons.reserve(100); // to avoid reallocation;
1185
1186     //set nodes shift within grid->_nodes from the node 000 
1187     size_t dx = _grid->NodeIndexDX();
1188     size_t dy = _grid->NodeIndexDY();
1189     size_t dz = _grid->NodeIndexDZ();
1190     size_t i000 = 0;
1191     size_t i100 = i000 + dx;
1192     size_t i010 = i000 + dy;
1193     size_t i110 = i010 + dx;
1194     size_t i001 = i000 + dz;
1195     size_t i101 = i100 + dz;
1196     size_t i011 = i010 + dz;
1197     size_t i111 = i110 + dz;
1198     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V000 )] = i000;
1199     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V100 )] = i100;
1200     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V010 )] = i010;
1201     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V110 )] = i110;
1202     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V001 )] = i001;
1203     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V101 )] = i101;
1204     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V011 )] = i011;
1205     _nodeShift[ SMESH_Block::ShapeIndex( SMESH_Block::ID_V111 )] = i111;
1206
1207     vector< int > idVec;
1208     // set nodes to links
1209     for ( int linkID = SMESH_Block::ID_Ex00; linkID <= SMESH_Block::ID_E11z; ++linkID )
1210     {
1211       SMESH_Block::GetEdgeVertexIDs( linkID, idVec );
1212       _Link& link = _hexLinks[ SMESH_Block::ShapeIndex( linkID )];
1213       link._nodes[0] = &_hexNodes[ SMESH_Block::ShapeIndex( idVec[0] )];
1214       link._nodes[1] = &_hexNodes[ SMESH_Block::ShapeIndex( idVec[1] )];
1215       link._intNodes.reserve( 10 ); // to avoid reallocation
1216       link._splits.reserve( 10 );
1217     }
1218
1219     // set links to faces
1220     int interlace[4] = { 0, 3, 1, 2 }; // to walk by links around a face: { u0, 1v, u1, 0v }
1221     for ( int faceID = SMESH_Block::ID_Fxy0; faceID <= SMESH_Block::ID_F1yz; ++faceID )
1222     {
1223       SMESH_Block::GetFaceEdgesIDs( faceID, idVec );
1224       _Face& quad = _hexQuads[ SMESH_Block::ShapeIndex( faceID )];
1225       bool revFace = ( faceID == SMESH_Block::ID_Fxy0 ||
1226                        faceID == SMESH_Block::ID_Fx1z ||
1227                        faceID == SMESH_Block::ID_F0yz );
1228       quad._links.resize(4);
1229       vector<_OrientedLink>::iterator         frwLinkIt = quad._links.begin();
1230       vector<_OrientedLink>::reverse_iterator revLinkIt = quad._links.rbegin();
1231       for ( int i = 0; i < 4; ++i )
1232       {
1233         bool revLink = revFace;
1234         if ( i > 1 ) // reverse links u1 and v0
1235           revLink = !revLink;
1236         _OrientedLink& link = revFace ? *revLinkIt++ : *frwLinkIt++;
1237         link = _OrientedLink( & _hexLinks[ SMESH_Block::ShapeIndex( idVec[interlace[i]] )],
1238                               revLink );
1239       }
1240     }
1241   }
1242   //================================================================================
1243   /*!
1244    * \brief Copy constructor
1245    */
1246   Hexahedron::Hexahedron( const Hexahedron& other )
1247     :_grid( other._grid ), _sizeThreshold( other._sizeThreshold ), _nbIntNodes(0)
1248   {
1249     _polygons.reserve(100); // to avoid reallocation;
1250
1251     for ( int i = 0; i < 8; ++i )
1252       _nodeShift[i] = other._nodeShift[i];
1253
1254     for ( int i = 0; i < 12; ++i )
1255     {
1256       const _Link& srcLink = other._hexLinks[ i ];
1257       _Link&       tgtLink = this->_hexLinks[ i ];
1258       tgtLink._nodes[0] = _hexNodes + ( srcLink._nodes[0] - other._hexNodes );
1259       tgtLink._nodes[1] = _hexNodes + ( srcLink._nodes[1] - other._hexNodes );
1260       tgtLink._intNodes.reserve( 10 ); // to avoid reallocation
1261       tgtLink._splits.reserve( 10 );
1262     }
1263
1264     for ( int i = 0; i < 6; ++i )
1265     {
1266       const _Face& srcQuad = other._hexQuads[ i ];
1267       _Face&       tgtQuad = this->_hexQuads[ i ];
1268       tgtQuad._links.resize(4);
1269       for ( int j = 0; j < 4; ++j )
1270       {
1271         const _OrientedLink& srcLink = srcQuad._links[ j ];
1272         _OrientedLink&       tgtLink = tgtQuad._links[ j ];
1273         tgtLink._reverse = srcLink._reverse;
1274         tgtLink._link    = _hexLinks + ( srcLink._link - other._hexLinks );
1275       }
1276     }
1277   }
1278   
1279   //================================================================================
1280   /*!
1281    * \brief Initializes its data by given grid cell
1282    */
1283   void Hexahedron::init( size_t i, size_t j, size_t k )
1284   {
1285     _i = i; _j = j; _k = k;
1286     // set nodes of grid to nodes of the hexahedron and
1287     // count nodes at hexahedron corners located IN and ON geometry
1288     _nbCornerNodes = _nbBndNodes = 0;
1289     _origNodeInd   = _grid->NodeIndex( i,j,k );
1290     for ( int iN = 0; iN < 8; ++iN )
1291     {
1292       _hexNodes[iN]._node = _grid->_nodes[ _origNodeInd + _nodeShift[iN] ];
1293       _nbCornerNodes += bool( _hexNodes[iN]._node );
1294       _nbBndNodes    += _grid->_isBndNode[ _origNodeInd + _nodeShift[iN] ];
1295     }
1296
1297     _sideLength[0] = _grid->_coords[0][i+1] - _grid->_coords[0][i];
1298     _sideLength[1] = _grid->_coords[1][j+1] - _grid->_coords[1][j];
1299     _sideLength[2] = _grid->_coords[2][k+1] - _grid->_coords[2][k];
1300
1301     if ( _nbCornerNodes < 8 && _nbIntNodes + _nbCornerNodes > 3)
1302     {
1303       _Link split;
1304       // create sub-links (_splits) by splitting links with _intNodes
1305       for ( int iLink = 0; iLink < 12; ++iLink )
1306       {
1307         _Link& link = _hexLinks[ iLink ];
1308         link._splits.clear();
1309         split._nodes[ 0 ] = link._nodes[0];
1310         for ( size_t i = 0; i < link._intNodes.size(); ++ i )
1311         {
1312           if ( split._nodes[ 0 ]->Node() )
1313           {
1314             split._nodes[ 1 ] = &link._intNodes[i];
1315             link._splits.push_back( split );
1316           }
1317           split._nodes[ 0 ] = &link._intNodes[i];
1318         }
1319         if ( link._nodes[ 1 ]->Node() && split._nodes[ 0 ]->Node() )
1320         {
1321           split._nodes[ 1 ] = link._nodes[1];
1322           link._splits.push_back( split );
1323         }
1324       }
1325     }
1326   }
1327   //================================================================================
1328   /*!
1329    * \brief Initializes its data by given grid cell (countered from zero)
1330    */
1331   void Hexahedron::init( size_t iCell )
1332   {
1333     size_t iNbCell = _grid->_coords[0].size() - 1;
1334     size_t jNbCell = _grid->_coords[1].size() - 1;
1335     _i = iCell % iNbCell;
1336     _j = ( iCell % ( iNbCell * jNbCell )) / iNbCell;
1337     _k = iCell / iNbCell / jNbCell;
1338     init( _i, _j, _k );
1339   }
1340
1341   //================================================================================
1342   /*!
1343    * \brief Compute mesh volumes resulted from intersection of the Hexahedron
1344    */
1345   void Hexahedron::ComputeElements()
1346   {
1347     Init();
1348
1349     if ( _nbCornerNodes + _nbIntNodes < 4 )
1350       return;
1351
1352     if ( _nbBndNodes == _nbCornerNodes && isInHole() )
1353       return;
1354
1355     _polygons.clear();
1356
1357     vector<const SMDS_MeshNode* > polyhedraNodes;
1358     vector<int>                   quantities;
1359
1360     // create polygons from quadrangles and get their nodes
1361
1362     vector<_Node*> nodes;
1363     nodes.reserve( _nbCornerNodes + _nbIntNodes );
1364
1365     _Link polyLink;
1366     polyLink._faces.reserve( 1 );
1367
1368     for ( int iF = 0; iF < 6; ++iF ) // loop on 6 sides of a hexahedron
1369     {
1370       const _Face& quad = _hexQuads[ iF ] ;
1371
1372       _polygons.resize( _polygons.size() + 1 );
1373       _Face& polygon = _polygons.back();
1374       polygon._links.clear();
1375       polygon._polyLinks.clear(); polygon._polyLinks.reserve( 10 );
1376
1377       // add splits of a link to a polygon and collect info on nodes
1378       //int nbIn = 0, nbOut = 0, nbCorners = 0;
1379       nodes.clear();
1380       for ( int iE = 0; iE < 4; ++iE ) // loop on 4 sides of a quadrangle
1381       {
1382         int nbSpits = quad._links[ iE ].NbResultLinks();
1383         for ( int iS = 0; iS < nbSpits; ++iS )
1384         {
1385           _OrientedLink split = quad._links[ iE ].ResultLink( iS );
1386           _Node* n = split.FirstNode();
1387           if ( !polygon._links.empty() )
1388           {
1389             _Node* nPrev = polygon._links.back().LastNode();
1390             if ( nPrev != n )
1391             {
1392               polyLink._nodes[0] = nPrev;
1393               polyLink._nodes[1] = n;
1394               polygon._polyLinks.push_back( polyLink );
1395               polygon._links.push_back( _OrientedLink( &polygon._polyLinks.back() ));
1396               nodes.push_back( nPrev );
1397             }
1398           }
1399           polygon._links.push_back( split );
1400           nodes.push_back( n );
1401         }
1402       }
1403       if ( polygon._links.size() > 1 )
1404       {
1405         _Node* n1 = polygon._links.back().LastNode();
1406         _Node* n2 = polygon._links.front().FirstNode();
1407         if ( n1 != n2 )
1408         {
1409           polyLink._nodes[0] = n1;
1410           polyLink._nodes[1] = n2;
1411           polygon._polyLinks.push_back( polyLink );
1412           polygon._links.push_back( _OrientedLink( &polygon._polyLinks.back() ));
1413           nodes.push_back( n1 );
1414         }
1415         // add polygon to its links
1416         for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
1417           polygon._links[ iL ]._link->_faces.push_back( &polygon );
1418         // store polygon nodes
1419         quantities.push_back( nodes.size() );
1420         for ( size_t i = 0; i < nodes.size(); ++i )
1421           polyhedraNodes.push_back( nodes[i]->Node() );
1422       }
1423       else
1424       {
1425         _polygons.resize( _polygons.size() - 1 );
1426       }
1427     }
1428
1429     // create polygons closing holes in a polyhedron
1430
1431     // find free links
1432     vector< _OrientedLink* > freeLinks;
1433     for ( size_t iP = 0; iP < _polygons.size(); ++iP )
1434     {
1435       _Face& polygon = _polygons[ iP ];
1436       for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
1437         if ( polygon._links[ iL ]._link->_faces.size() < 2 )
1438           freeLinks.push_back( & polygon._links[ iL ]);
1439     }
1440     // make closed chains of free links
1441     int nbFreeLinks = freeLinks.size();
1442     if ( 0 < nbFreeLinks && nbFreeLinks < 3 ) return;
1443     while ( nbFreeLinks > 0 )
1444     {
1445       nodes.clear();
1446       _polygons.resize( _polygons.size() + 1 );
1447       _Face& polygon = _polygons.back();
1448       polygon._links.clear();
1449
1450       // get a remaining link to start from
1451       _OrientedLink* curLink = 0;
1452       for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
1453         if (( curLink = freeLinks[ iL ] ))
1454           freeLinks[ iL ] = 0;
1455       nodes.push_back( curLink->LastNode() );
1456       polygon._links.push_back( *curLink );
1457
1458       // find all links connected to curLink
1459       _Node* curNode = 0;
1460       do
1461       {
1462         curNode = curLink->FirstNode();
1463         curLink = 0;
1464         for ( size_t iL = 0; iL < freeLinks.size() && !curLink; ++iL )
1465           if ( freeLinks[ iL ] && freeLinks[ iL ]->LastNode() == curNode )
1466           {
1467             curLink = freeLinks[ iL ];
1468             freeLinks[ iL ] = 0;
1469             nodes.push_back( curNode );
1470             polygon._links.push_back( *curLink );
1471           }
1472       } while ( curLink );
1473
1474       nbFreeLinks -= polygon._links.size();
1475
1476       if ( curNode != nodes.front() || polygon._links.size() < 3 )
1477         return; // closed polygon not found -> invalid polyhedron
1478
1479       quantities.push_back( nodes.size() );
1480       for ( size_t i = 0; i < nodes.size(); ++i )
1481         polyhedraNodes.push_back( nodes[i]->Node() );
1482
1483       // add polygon to its links and reverse links
1484       for ( size_t i = 0; i < polygon._links.size(); ++i )
1485       {
1486         polygon._links[i].Reverse();
1487         polygon._links[i]._link->_faces.push_back( &polygon );
1488       }
1489
1490       //const size_t firstPoly = _polygons.size();
1491     }
1492
1493     if ( ! checkPolyhedronSize() )
1494     {
1495       return;
1496     }
1497
1498     // create a classic cell if possible
1499     const int nbNodes = _nbCornerNodes + _nbIntNodes;
1500     bool isClassicElem = false;
1501     if (      nbNodes == 8 && _polygons.size() == 6 ) isClassicElem = addHexa();
1502     else if ( nbNodes == 4 && _polygons.size() == 4 ) isClassicElem = addTetra();
1503     else if ( nbNodes == 6 && _polygons.size() == 5 ) isClassicElem = addPenta();
1504     else if ( nbNodes == 5 && _polygons.size() == 5 ) isClassicElem = addPyra ();
1505     if ( !isClassicElem )
1506       _volumeDefs.set( polyhedraNodes, quantities );
1507   }
1508   //================================================================================
1509   /*!
1510    * \brief Create elements in the mesh
1511    */
1512   int Hexahedron::MakeElements(SMESH_MesherHelper& helper)
1513   {
1514     SMESHDS_Mesh* mesh = helper.GetMeshDS();
1515
1516     size_t nbCells[3] = { _grid->_coords[0].size() - 1,
1517                           _grid->_coords[1].size() - 1,
1518                           _grid->_coords[2].size() - 1 };
1519     const size_t nbGridCells = nbCells[0] *nbCells [1] * nbCells[2];
1520     vector< Hexahedron* > intersectedHex( nbGridCells, 0 );
1521     int nbIntHex = 0;
1522
1523     // set intersection nodes from GridLine's to links of intersectedHex
1524     int i,j,k, iDirOther[3][2] = {{ 1,2 },{ 0,2 },{ 0,1 }};
1525     for ( int iDir = 0; iDir < 3; ++iDir )
1526     {
1527       int dInd[4][3] = { {0,0,0}, {0,0,0}, {0,0,0}, {0,0,0} };
1528       dInd[1][ iDirOther[iDir][0] ] = -1;
1529       dInd[2][ iDirOther[iDir][1] ] = -1;
1530       dInd[3][ iDirOther[iDir][0] ] = -1; dInd[3][ iDirOther[iDir][1] ] = -1;
1531       // loop on GridLine's parallel to iDir
1532       LineIndexer lineInd = _grid->GetLineIndexer( iDir );
1533       for ( ; lineInd.More(); ++lineInd )
1534       {
1535         GridLine& line = _grid->_lines[ iDir ][ lineInd.LineIndex() ];
1536         multiset< IntersectionPoint >::const_iterator ip = line._intPoints.begin();
1537         for ( ; ip != line._intPoints.end(); ++ip )
1538         {
1539           if ( !ip->_node ) continue;
1540           lineInd.SetIndexOnLine( ip->_indexOnLine );
1541           for ( int iL = 0; iL < 4; ++iL ) // loop on 4 cells sharing a link
1542           {
1543             i = int(lineInd.I()) + dInd[iL][0];
1544             j = int(lineInd.J()) + dInd[iL][1];
1545             k = int(lineInd.K()) + dInd[iL][2];
1546             if ( i < 0 || i >= nbCells[0] ||
1547                  j < 0 || j >= nbCells[1] ||
1548                  k < 0 || k >= nbCells[2] ) continue;
1549
1550             const size_t hexIndex = _grid->CellIndex( i,j,k );
1551             Hexahedron *& hex = intersectedHex[ hexIndex ];
1552             if ( !hex)
1553             {
1554               hex = new Hexahedron( *this );
1555               hex->_i = i;
1556               hex->_j = j;
1557               hex->_k = k;
1558               ++nbIntHex;
1559             }
1560             const int iLink = iL + iDir * 4;
1561             hex->_hexLinks[iLink]._intNodes.push_back( _Node( 0, &(*ip) ));
1562             hex->_nbIntNodes++;
1563           }
1564         }
1565       }
1566     }
1567
1568     // add not split hexadrons to the mesh
1569     int nbAdded = 0;
1570     vector<int> intHexInd( nbIntHex );
1571     nbIntHex = 0;
1572     for ( size_t i = 0; i < intersectedHex.size(); ++i )
1573     {
1574       Hexahedron * & hex = intersectedHex[ i ];
1575       if ( hex )
1576       {
1577         intHexInd[ nbIntHex++ ] = i;
1578         if ( hex->_nbIntNodes > 0 ) continue;
1579         init( hex->_i, hex->_j, hex->_k );
1580       }
1581       else
1582       {    
1583         init( i );
1584       }
1585       if ( _nbCornerNodes == 8 && ( _nbBndNodes < _nbCornerNodes || !isInHole() ))
1586       {
1587         // order of _hexNodes is defined by enum SMESH_Block::TShapeID
1588         SMDS_MeshElement* el =
1589           mesh->AddVolume( _hexNodes[0].Node(), _hexNodes[2].Node(),
1590                            _hexNodes[3].Node(), _hexNodes[1].Node(),
1591                            _hexNodes[4].Node(), _hexNodes[6].Node(),
1592                            _hexNodes[7].Node(), _hexNodes[5].Node() );
1593         mesh->SetMeshElementOnShape( el, helper.GetSubShapeID() );
1594         ++nbAdded;
1595         if ( hex )
1596         {
1597           delete hex;
1598           intersectedHex[ i ] = 0;
1599           --nbIntHex;
1600         }
1601       }
1602       else if ( _nbCornerNodes > 3  && !hex )
1603       {
1604         // all intersection of hex with geometry are at grid nodes
1605         hex = new Hexahedron( *this );
1606         hex->init( i );
1607         intHexInd.push_back(0);
1608         intHexInd[ nbIntHex++ ] = i;
1609       }
1610     }
1611
1612     // add elements resulted from hexadron intersection
1613 #ifdef WITH_TBB
1614     intHexInd.resize( nbIntHex );
1615     tbb::parallel_for ( tbb::blocked_range<size_t>( 0, nbIntHex ),
1616                         ParallelHexahedron( intersectedHex, intHexInd ),
1617                         tbb::simple_partitioner()); // ComputeElements() is called here
1618     for ( size_t i = 0; i < intHexInd.size(); ++i )
1619       if ( Hexahedron * hex = intersectedHex[ intHexInd[ i ]] )
1620         nbAdded += hex->addElements( helper );
1621 #else
1622     for ( size_t i = 0; i < intHexInd.size(); ++i )
1623       if ( Hexahedron * hex = intersectedHex[ intHexInd[ i ]] )
1624       {
1625         hex->ComputeElements();
1626         nbAdded += hex->addElements( helper );
1627       }
1628 #endif
1629
1630     for ( size_t i = 0; i < intersectedHex.size(); ++i )
1631       if ( intersectedHex[ i ] )
1632         delete intersectedHex[ i ];
1633
1634     return nbAdded;
1635   }
1636
1637   //================================================================================
1638   /*!
1639    * \brief Adds computed elements to the mesh
1640    */
1641   int Hexahedron::addElements(SMESH_MesherHelper& helper)
1642   {
1643     int nbAdded = 0;
1644     // add elements resulted from hexahedron intersection
1645     //for ( size_t i = 0; i < _volumeDefs.size(); ++i )
1646     {
1647       vector< const SMDS_MeshNode* >& nodes = _volumeDefs._nodes;
1648       
1649       if ( !_volumeDefs._quantities.empty() )
1650       {
1651         helper.AddPolyhedralVolume( nodes, _volumeDefs._quantities );
1652       }
1653       else
1654       {
1655         switch ( nodes.size() )
1656         {
1657         case 8: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],
1658                                   nodes[4],nodes[5],nodes[6],nodes[7] );
1659           break;
1660         case 4: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3] );
1661           break;
1662         case 6: helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3], nodes[4],nodes[5] );
1663           break;
1664         case 5:
1665           helper.AddVolume( nodes[0],nodes[1],nodes[2],nodes[3],nodes[4] );
1666           break;
1667         }
1668       }
1669       nbAdded += int ( _volumeDefs._nodes.size() > 0 );
1670     }
1671
1672     return nbAdded;
1673   }
1674   //================================================================================
1675   /*!
1676    * \brief Return true if the element is in a hole
1677    */
1678   bool Hexahedron::isInHole() const
1679   {
1680     const int ijk[3] = { _i, _j, _k };
1681     IntersectionPoint curIntPnt;
1682
1683     // consider a cell to be in a hole if all links in any direction
1684     // comes OUT of geometry
1685     for ( int iDir = 0; iDir < 3; ++iDir )
1686     {
1687       const vector<double>& coords = _grid->_coords[ iDir ];
1688       LineIndexer               li = _grid->GetLineIndexer( iDir );
1689       li.SetIJK( _i,_j,_k );
1690       size_t lineIndex[4] = { li.LineIndex  (),
1691                               li.LineIndex10(),
1692                               li.LineIndex01(),
1693                               li.LineIndex11() };
1694       bool allLinksOut = true, hasLinks = false;
1695       for ( int iL = 0; iL < 4 && allLinksOut; ++iL ) // loop on 4 links parallel to iDir
1696       {
1697         const _Link& link = _hexLinks[ iL + 4*iDir ];
1698         // check transition of the first node of a link
1699         const IntersectionPoint* firstIntPnt = 0;
1700         if ( link._nodes[0]->Node() ) // 1st node is a hexa corner
1701         {
1702           curIntPnt._paramOnLine = coords[ ijk[ iDir ]] - coords[0];
1703           const GridLine& line = _grid->_lines[ iDir ][ lineIndex[ iL ]];
1704           multiset< IntersectionPoint >::const_iterator ip =
1705             line._intPoints.upper_bound( curIntPnt );
1706           --ip;
1707           firstIntPnt = &(*ip);
1708         }
1709         else if ( !link._intNodes.empty() )
1710         {
1711           firstIntPnt = link._intNodes[0]._intPoint;
1712         }
1713
1714         if ( firstIntPnt )
1715         {
1716           hasLinks = true;
1717           allLinksOut = ( firstIntPnt->_transition == Trans_OUT );
1718         }
1719       }
1720       if ( hasLinks && allLinksOut )
1721         return true;
1722     }
1723     return false;
1724   }
1725
1726   //================================================================================
1727   /*!
1728    * \brief Return true if a polyhedron passes _sizeThreshold criterion
1729    */
1730   bool Hexahedron::checkPolyhedronSize() const
1731   {
1732     double volume = 0;
1733     for ( size_t iP = 0; iP < _polygons.size(); ++iP )
1734     {
1735       const _Face& polygon = _polygons[iP];
1736       gp_XYZ area (0,0,0);
1737       SMESH_TNodeXYZ p1 ( polygon._links[ 0 ].FirstNode()->Node() );
1738       for ( size_t iL = 0; iL < polygon._links.size(); ++iL )
1739       {
1740         SMESH_TNodeXYZ p2 ( polygon._links[ iL ].LastNode()->Node() );
1741         area += p1 ^ p2;
1742         p1 = p2;
1743       }
1744       volume += p1 * area;
1745     }
1746     volume /= 6;
1747
1748     double initVolume = _sideLength[0] * _sideLength[1] * _sideLength[2];
1749
1750     return volume > initVolume / _sizeThreshold;
1751   }
1752   //================================================================================
1753   /*!
1754    * \brief Tries to create a hexahedron
1755    */
1756   bool Hexahedron::addHexa()
1757   {
1758     if ( _polygons[0]._links.size() != 4 ||
1759          _polygons[1]._links.size() != 4 ||
1760          _polygons[2]._links.size() != 4 ||
1761          _polygons[3]._links.size() != 4 ||
1762          _polygons[4]._links.size() != 4 ||
1763          _polygons[5]._links.size() != 4   )
1764       return false;
1765     const SMDS_MeshNode* nodes[8];
1766     int nbN = 0;
1767     for ( int iL = 0; iL < 4; ++iL )
1768     {
1769       // a base node
1770       nodes[iL] = _polygons[0]._links[iL].FirstNode()->Node();
1771       ++nbN;
1772
1773       // find a top node above the base node
1774       _Link* link = _polygons[0]._links[iL]._link;
1775       ASSERT( link->_faces.size() > 1 );
1776       // a quadrangle sharing <link> with _polygons[0]
1777       _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[0] )];
1778       for ( int i = 0; i < 4; ++i )
1779         if ( quad->_links[i]._link == link )
1780         {
1781           // 1st node of a link opposite to <link> in <quad>
1782           nodes[iL+4] = quad->_links[(i+2)%4].FirstNode()->Node();
1783           ++nbN;
1784           break;
1785         }
1786     }
1787     if ( nbN == 8 )
1788       _volumeDefs.set( vector< const SMDS_MeshNode* >( nodes, nodes+8 ));
1789
1790     return nbN == 8;
1791   }
1792   //================================================================================
1793   /*!
1794    * \brief Tries to create a tetrahedron
1795    */
1796   bool Hexahedron::addTetra()
1797   {
1798     const SMDS_MeshNode* nodes[4];
1799     nodes[0] = _polygons[0]._links[0].FirstNode()->Node();
1800     nodes[1] = _polygons[0]._links[1].FirstNode()->Node();
1801     nodes[2] = _polygons[0]._links[2].FirstNode()->Node();
1802
1803     _Link* link = _polygons[0]._links[0]._link;
1804     ASSERT( link->_faces.size() > 1 );
1805
1806     // a triangle sharing <link> with _polygons[0]
1807     _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[0] )];
1808     for ( int i = 0; i < 3; ++i )
1809       if ( tria->_links[i]._link == link )
1810       {
1811         nodes[3] = tria->_links[(i+1)%3].LastNode()->Node();
1812         _volumeDefs.set( vector< const SMDS_MeshNode* >( nodes, nodes+4 ));
1813         return true;
1814       }
1815
1816     return false;
1817   }
1818   //================================================================================
1819   /*!
1820    * \brief Tries to create a pentahedron
1821    */
1822   bool Hexahedron::addPenta()
1823   {
1824     // find a base triangular face
1825     int iTri = -1;
1826     for ( int iF = 0; iF < 5 && iTri < 0; ++iF )
1827       if ( _polygons[ iF ]._links.size() == 3 )
1828         iTri = iF;
1829     if ( iTri < 0 ) return false;
1830
1831     // find nodes
1832     const SMDS_MeshNode* nodes[6];
1833     int nbN = 0;
1834     for ( int iL = 0; iL < 3; ++iL )
1835     {
1836       // a base node
1837       nodes[iL] = _polygons[ iTri ]._links[iL].FirstNode()->Node();
1838       ++nbN;
1839
1840       // find a top node above the base node
1841       _Link* link = _polygons[ iTri ]._links[iL]._link;
1842       ASSERT( link->_faces.size() > 1 );
1843       // a quadrangle sharing <link> with a base triangle
1844       _Face* quad = link->_faces[ bool( link->_faces[0] == & _polygons[ iTri ] )];
1845       if ( quad->_links.size() != 4 ) return false;
1846       for ( int i = 0; i < 4; ++i )
1847         if ( quad->_links[i]._link == link )
1848         {
1849           // 1st node of a link opposite to <link> in <quad>
1850           nodes[iL+3] = quad->_links[(i+2)%4].FirstNode()->Node();
1851           ++nbN;
1852           break;
1853         }
1854     }
1855     if ( nbN == 6 )
1856       _volumeDefs.set( vector< const SMDS_MeshNode* >( nodes, nodes+6 ));
1857
1858     return ( nbN == 6 );
1859   }
1860   //================================================================================
1861   /*!
1862    * \brief Tries to create a pyramid
1863    */
1864   bool Hexahedron::addPyra()
1865   {
1866     // find a base quadrangle
1867     int iQuad = -1;
1868     for ( int iF = 0; iF < 5 && iQuad < 0; ++iF )
1869       if ( _polygons[ iF ]._links.size() == 4 )
1870         iQuad = iF;
1871     if ( iQuad < 0 ) return false;
1872
1873     // find nodes
1874     const SMDS_MeshNode* nodes[5];
1875     nodes[0] = _polygons[iQuad]._links[0].FirstNode()->Node();
1876     nodes[1] = _polygons[iQuad]._links[1].FirstNode()->Node();
1877     nodes[2] = _polygons[iQuad]._links[2].FirstNode()->Node();
1878     nodes[3] = _polygons[iQuad]._links[3].FirstNode()->Node();
1879
1880     _Link* link = _polygons[iQuad]._links[0]._link;
1881     ASSERT( link->_faces.size() > 1 );
1882
1883     // a triangle sharing <link> with a base quadrangle
1884     _Face* tria = link->_faces[ bool( link->_faces[0] == & _polygons[ iQuad ] )];
1885     if ( tria->_links.size() != 3 ) return false;
1886     for ( int i = 0; i < 3; ++i )
1887       if ( tria->_links[i]._link == link )
1888       {
1889         nodes[4] = tria->_links[(i+1)%3].LastNode()->Node();
1890         _volumeDefs.set( vector< const SMDS_MeshNode* >( nodes, nodes+5 ));
1891         return true;
1892       }
1893
1894     return false;
1895   }
1896
1897 } // namespace
1898
1899 //=============================================================================
1900 /*!
1901  * \brief Generates 3D structured Cartesian mesh in the internal part of
1902  * solid shapes and polyhedral volumes near the shape boundary.
1903  *  \param theMesh - mesh to fill in
1904  *  \param theShape - a compound of all SOLIDs to mesh
1905  *  \retval bool - true in case of success
1906  */
1907 //=============================================================================
1908
1909 bool StdMeshers_Cartesian_3D::Compute(SMESH_Mesh &         theMesh,
1910                                       const TopoDS_Shape & theShape)
1911 {
1912   // The algorithm generates the mesh in following steps:
1913
1914   // 1) Intersection of grid lines with the geometry boundary.
1915   // This step allows to find out if a given node of the initial grid is
1916   // inside or outside the geometry.
1917
1918   // 2) For each cell of the grid, check how many of it's nodes are outside
1919   // of the geometry boundary. Depending on a result of this check
1920   // - skip a cell, if all it's nodes are outside
1921   // - skip a cell, if it is too small according to the size threshold
1922   // - add a hexahedron in the mesh, if all nodes are inside
1923   // - add a polyhedron in the mesh, if some nodes are inside and some outside
1924
1925   _computeCanceled = false;
1926
1927   try
1928   {
1929     Grid grid;
1930
1931     TopTools_MapOfShape faceMap;
1932     for ( TopExp_Explorer fExp( theShape, TopAbs_FACE ); fExp.More(); fExp.Next() )
1933       if ( !faceMap.Add( fExp.Current() ))
1934         faceMap.Remove( fExp.Current() ); // remove a face shared by two solids
1935
1936     Bnd_Box shapeBox;
1937     vector<FaceGridIntersector> facesItersectors( faceMap.Extent() );
1938     TopTools_MapIteratorOfMapOfShape faceMppIt( faceMap );
1939     for ( int i = 0; faceMppIt.More(); faceMppIt.Next(), ++i )
1940     {
1941       facesItersectors[i]._face = TopoDS::Face( faceMppIt.Key() );
1942       facesItersectors[i]._grid = &grid;
1943       shapeBox.Add( facesItersectors[i].GetFaceBndBox() );
1944     }
1945
1946     vector<double> xCoords, yCoords, zCoords;
1947     _hyp->GetCoordinates( xCoords, yCoords, zCoords, shapeBox );
1948
1949     grid.SetCoordinates( xCoords, yCoords, zCoords, theShape );
1950
1951     // check if the grid encloses the shape
1952     if ( !_hyp->IsGridBySpacing(0) ||
1953          !_hyp->IsGridBySpacing(1) ||
1954          !_hyp->IsGridBySpacing(2) )
1955     {
1956       Bnd_Box gridBox;
1957       gridBox.Add( gp_Pnt( xCoords[0], yCoords[0], zCoords[0] ));
1958       gridBox.Add( gp_Pnt( xCoords.back(), yCoords.back(), zCoords.back() ));
1959       double x0,y0,z0, x1,y1,z1;
1960       shapeBox.Get(x0,y0,z0, x1,y1,z1);
1961       if ( gridBox.IsOut( gp_Pnt( x0,y0,z0 )) ||
1962            gridBox.IsOut( gp_Pnt( x1,y1,z1 )))
1963         for ( size_t i = 0; i < facesItersectors.size(); ++i )
1964         {
1965           if ( !facesItersectors[i].IsInGrid( gridBox ))
1966             return error("The grid doesn't enclose the geometry");
1967 #ifdef ELLIPSOLID_WORKAROUND
1968           delete facesItersectors[i]._surfaceInt, facesItersectors[i]._surfaceInt = 0;
1969 #endif
1970         }
1971     }
1972     if ( _computeCanceled ) return false;
1973
1974 #ifdef WITH_TBB
1975     { // copy partner faces and curves of not thread-safe types
1976       set< const Standard_Transient* > tshapes;
1977       BRepBuilderAPI_Copy copier;
1978       for ( size_t i = 0; i < facesItersectors.size(); ++i )
1979       {
1980         if ( !facesItersectors[i].IsThreadSafe(tshapes) )
1981         {
1982           copier.Perform( facesItersectors[i]._face );
1983           facesItersectors[i]._face = TopoDS::Face( copier );
1984         }
1985       }
1986     }
1987     // Intersection of grid lines with the geometry boundary.
1988     tbb::parallel_for ( tbb::blocked_range<size_t>( 0, facesItersectors.size() ),
1989                         ParallelIntersector( facesItersectors ),
1990                         tbb::simple_partitioner());
1991 #else
1992     for ( size_t i = 0; i < facesItersectors.size(); ++i )
1993       facesItersectors[i].Intersect();
1994 #endif
1995
1996     // put interesection points onto the GridLine's; this is done after intersection
1997     // to avoid contention of facesItersectors for writing into the same GridLine
1998     // in case of parallel work of facesItersectors
1999     for ( size_t i = 0; i < facesItersectors.size(); ++i )
2000       facesItersectors[i].StoreIntersections();
2001
2002     SMESH_MesherHelper helper( theMesh );
2003     TopExp_Explorer solidExp (theShape, TopAbs_SOLID);
2004     helper.SetSubShape( solidExp.Current() );
2005     helper.SetElementsOnShape( true );
2006
2007     if ( _computeCanceled ) return false;
2008
2009     // create nodes on the geometry
2010     grid.ComputeNodes(helper);
2011
2012     if ( _computeCanceled ) return false;
2013
2014     // create volume elements
2015     Hexahedron hex( _hyp->GetSizeThreshold(), &grid );
2016     int nbAdded = hex.MakeElements( helper );
2017
2018     SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
2019     if ( nbAdded > 0 )
2020     {
2021       // make all SOLIDS computed
2022       if ( SMESHDS_SubMesh* sm1 = meshDS->MeshElements( solidExp.Current()) )
2023       {
2024         SMDS_ElemIteratorPtr volIt = sm1->GetElements();
2025         for ( ; solidExp.More() && volIt->more(); solidExp.Next() )
2026         {
2027           const SMDS_MeshElement* vol = volIt->next();
2028           sm1->RemoveElement( vol, /*isElemDeleted=*/false );
2029           meshDS->SetMeshElementOnShape( vol, solidExp.Current() );
2030         }
2031       }
2032       // make other sub-shapes computed
2033       setSubmeshesComputed( theMesh, theShape );
2034     }
2035
2036     // remove free nodes
2037     if ( SMESHDS_SubMesh * smDS = meshDS->MeshElements( helper.GetSubShapeID() ))
2038     {
2039       // intersection nodes
2040       for ( int iDir = 0; iDir < 3; ++iDir )
2041       {
2042         vector< GridLine >& lines = grid._lines[ iDir ];
2043         for ( size_t i = 0; i < lines.size(); ++i )
2044         {
2045           multiset< IntersectionPoint >::iterator ip = lines[i]._intPoints.begin();
2046           for ( ; ip != lines[i]._intPoints.end(); ++ip )
2047             if ( ip->_node && ip->_node->NbInverseElements() == 0 )
2048               meshDS->RemoveFreeNode( ip->_node, smDS, /*fromGroups=*/false );
2049         }
2050       }
2051       // grid nodes
2052       for ( size_t i = 0; i < grid._nodes.size(); ++i )
2053         if ( !grid._isBndNode[i] ) // nodes on boundary are already removed
2054           if ( grid._nodes[i] && grid._nodes[i]->NbInverseElements() == 0 )
2055             meshDS->RemoveFreeNode( grid._nodes[i], smDS, /*fromGroups=*/false );
2056     }
2057
2058     return nbAdded;
2059
2060   }
2061   // SMESH_ComputeError is not caught at SMESH_submesh level for an unknown reason
2062   catch ( SMESH_ComputeError& e)
2063   {
2064     return error( SMESH_ComputeErrorPtr( new SMESH_ComputeError( e )));
2065   }
2066   return false;
2067 }
2068
2069 //=============================================================================
2070 /*!
2071  *  Evaluate
2072  */
2073 //=============================================================================
2074
2075 bool StdMeshers_Cartesian_3D::Evaluate(SMESH_Mesh &         theMesh,
2076                                        const TopoDS_Shape & theShape,
2077                                        MapShapeNbElems&     theResMap)
2078 {
2079   // TODO
2080 //   std::vector<int> aResVec(SMDSEntity_Last);
2081 //   for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
2082 //   if(IsQuadratic) {
2083 //     aResVec[SMDSEntity_Quad_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
2084 //     int nb1d_face0_int = ( nb2d_face0*4 - nb1d ) / 2;
2085 //     aResVec[SMDSEntity_Node] = nb0d_face0 * ( 2*nb2d/nb1d - 1 ) - nb1d_face0_int * nb2d/nb1d;
2086 //   }
2087 //   else {
2088 //     aResVec[SMDSEntity_Node] = nb0d_face0 * ( nb2d/nb1d - 1 );
2089 //     aResVec[SMDSEntity_Cartesian] = nb2d_face0 * ( nb2d/nb1d );
2090 //   }
2091 //   SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
2092 //   aResMap.insert(std::make_pair(sm,aResVec));
2093
2094   return true;
2095 }
2096
2097 //=============================================================================
2098 namespace
2099 {
2100   /*!
2101    * \brief Event listener setting/unsetting _alwaysComputed flag to
2102    *        submeshes of inferior levels to prevent their computing
2103    */
2104   struct _EventListener : public SMESH_subMeshEventListener
2105   {
2106     string _algoName;
2107
2108     _EventListener(const string& algoName):
2109       SMESH_subMeshEventListener(/*isDeletable=*/true,"StdMeshers_Cartesian_3D::_EventListener"),
2110       _algoName(algoName)
2111     {}
2112     // --------------------------------------------------------------------------------
2113     // setting/unsetting _alwaysComputed flag to submeshes of inferior levels
2114     //
2115     static void setAlwaysComputed( const bool     isComputed,
2116                                    SMESH_subMesh* subMeshOfSolid)
2117     {
2118       SMESH_subMeshIteratorPtr smIt =
2119         subMeshOfSolid->getDependsOnIterator(/*includeSelf=*/false, /*complexShapeFirst=*/false);
2120       while ( smIt->more() )
2121       {
2122         SMESH_subMesh* sm = smIt->next();
2123         sm->SetIsAlwaysComputed( isComputed );
2124       }
2125     }
2126
2127     // --------------------------------------------------------------------------------
2128     // unsetting _alwaysComputed flag if "Cartesian_3D" was removed
2129     //
2130     virtual void ProcessEvent(const int          event,
2131                               const int          eventType,
2132                               SMESH_subMesh*     subMeshOfSolid,
2133                               SMESH_subMeshEventListenerData* data,
2134                               const SMESH_Hypothesis*         hyp = 0)
2135     {
2136       if ( eventType == SMESH_subMesh::COMPUTE_EVENT )
2137       {
2138         setAlwaysComputed( subMeshOfSolid->GetComputeState() == SMESH_subMesh::COMPUTE_OK,
2139                            subMeshOfSolid );
2140       }
2141       else
2142       {
2143         SMESH_Algo* algo3D = subMeshOfSolid->GetAlgo();
2144         if ( !algo3D || _algoName != algo3D->GetName() )
2145           setAlwaysComputed( false, subMeshOfSolid );
2146       }
2147     }
2148
2149     // --------------------------------------------------------------------------------
2150     // set the event listener
2151     //
2152     static void SetOn( SMESH_subMesh* subMeshOfSolid, const string& algoName )
2153     {
2154       subMeshOfSolid->SetEventListener( new _EventListener( algoName ),
2155                                         /*data=*/0,
2156                                         subMeshOfSolid );
2157     }
2158
2159   }; // struct _EventListener
2160
2161 } // namespace
2162
2163 //================================================================================
2164 /*!
2165  * \brief Sets event listener to submeshes if necessary
2166  *  \param subMesh - submesh where algo is set
2167  * This method is called when a submesh gets HYP_OK algo_state.
2168  * After being set, event listener is notified on each event of a submesh.
2169  */
2170 //================================================================================
2171
2172 void StdMeshers_Cartesian_3D::SetEventListener(SMESH_subMesh* subMesh)
2173 {
2174   _EventListener::SetOn( subMesh, GetName() );
2175 }
2176
2177 //================================================================================
2178 /*!
2179  * \brief Set _alwaysComputed flag to submeshes of inferior levels to avoid their computing
2180  */
2181 //================================================================================
2182
2183 void StdMeshers_Cartesian_3D::setSubmeshesComputed(SMESH_Mesh&         theMesh,
2184                                                    const TopoDS_Shape& theShape)
2185 {
2186   for ( TopExp_Explorer soExp( theShape, TopAbs_SOLID ); soExp.More(); soExp.Next() )
2187     _EventListener::setAlwaysComputed( true, theMesh.GetSubMesh( soExp.Current() ));
2188 }
2189