Salome HOME
2fbab5f95cb38fd2e781db67fda89752bcf15a6a
[modules/smesh.git] / src / StdMeshers / StdMeshers_ViscousLayers2D.cxx
1 // Copyright (C) 2007-2012  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 // File      : StdMeshers_ViscousLayers2D.cxx
21 // Created   : 23 Jul 2012
22 // Author    : Edward AGAPOV (eap)
23
24 #include "StdMeshers_ViscousLayers2D.hxx"
25
26 #include "SMDS_EdgePosition.hxx"
27 #include "SMDS_FaceOfNodes.hxx"
28 #include "SMDS_FacePosition.hxx"
29 #include "SMDS_MeshNode.hxx"
30 #include "SMDS_SetIterator.hxx"
31 #include "SMESHDS_Group.hxx"
32 #include "SMESHDS_Hypothesis.hxx"
33 #include "SMESH_Algo.hxx"
34 #include "SMESH_ComputeError.hxx"
35 #include "SMESH_ControlsDef.hxx"
36 #include "SMESH_Gen.hxx"
37 #include "SMESH_Group.hxx"
38 #include "SMESH_HypoFilter.hxx"
39 #include "SMESH_Mesh.hxx"
40 #include "SMESH_MesherHelper.hxx"
41 #include "SMESH_ProxyMesh.hxx"
42 #include "SMESH_Quadtree.hxx"
43 #include "SMESH_subMesh.hxx"
44 #include "SMESH_subMeshEventListener.hxx"
45 #include "StdMeshers_FaceSide.hxx"
46
47 #include "utilities.h"
48
49 #include <BRepAdaptor_Curve.hxx>
50 #include <BRepAdaptor_Curve2d.hxx>
51 #include <BRep_Tool.hxx>
52 #include <Bnd_B2d.hxx>
53 #include <Bnd_B3d.hxx>
54 #include <ElCLib.hxx>
55 #include <GCPnts_AbscissaPoint.hxx>
56 #include <Geom2dAdaptor_Curve.hxx>
57 #include <Geom2dInt_GInter.hxx>
58 #include <Geom2d_Circle.hxx>
59 #include <Geom2d_Line.hxx>
60 #include <Geom2d_TrimmedCurve.hxx>
61 #include <GeomAdaptor_Curve.hxx>
62 #include <Geom_Circle.hxx>
63 #include <Geom_Curve.hxx>
64 #include <Geom_Line.hxx>
65 #include <Geom_TrimmedCurve.hxx>
66 #include <IntRes2d_IntersectionPoint.hxx>
67 #include <Precision.hxx>
68 #include <Standard_ErrorHandler.hxx>
69 #include <TColStd_Array1OfReal.hxx>
70 #include <TopExp.hxx>
71 #include <TopExp_Explorer.hxx>
72 #include <TopTools_IndexedMapOfShape.hxx>
73 #include <TopTools_MapOfShape.hxx>
74 #include <TopoDS.hxx>
75 #include <TopoDS_Edge.hxx>
76 #include <TopoDS_Face.hxx>
77 #include <TopoDS_Vertex.hxx>
78 #include <gp_Ax1.hxx>
79 #include <gp_Vec.hxx>
80 #include <gp_XY.hxx>
81
82 #include <list>
83 #include <string>
84 #include <cmath>
85 #include <limits>
86
87 #define __myDEBUG
88
89 using namespace std;
90
91 //================================================================================
92 namespace VISCOUS_2D
93 {
94   typedef int TGeomID;
95
96   //--------------------------------------------------------------------------------
97   /*!
98    * \brief Proxy Mesh of FACE with viscous layers. It's needed only to 
99    *        redefine newSubmesh().
100    */
101   struct _ProxyMeshOfFace : public SMESH_ProxyMesh
102   {
103     //---------------------------------------------------
104     // Proxy sub-mesh of an EDGE. It contains nodes in _uvPtStructVec.
105     struct _EdgeSubMesh : public SMESH_ProxyMesh::SubMesh
106     {
107       _EdgeSubMesh(int index=0): SubMesh(index) {}
108       //virtual int NbElements() const { return _elements.size()+1; }
109       virtual int NbNodes() const { return Max( 0, _uvPtStructVec.size()-2 ); }
110       void SetUVPtStructVec(UVPtStructVec& vec) { _uvPtStructVec.swap( vec ); }
111     };
112     _ProxyMeshOfFace(const SMESH_Mesh& mesh): SMESH_ProxyMesh(mesh) {}
113     _EdgeSubMesh* GetEdgeSubMesh(int ID) { return (_EdgeSubMesh*) getProxySubMesh(ID); }
114     virtual SubMesh* newSubmesh(int index=0) const { return new _EdgeSubMesh(index); }
115   };
116   //--------------------------------------------------------------------------------
117   /*!
118    * \brief SMESH_subMeshEventListener used to store _ProxyMeshOfFace, computed
119    *        by _ViscousBuilder2D, in a SMESH_subMesh of the FACE.
120    *        This is to delete _ProxyMeshOfFace when StdMeshers_ViscousLayers2D
121    *        hypothesis is modified
122    */
123   struct _ProxyMeshHolder : public SMESH_subMeshEventListener
124   {
125     _ProxyMeshHolder( const TopoDS_Face&    face,
126                       SMESH_ProxyMesh::Ptr& mesh)
127       : SMESH_subMeshEventListener( /*deletable=*/true, Name() )
128     {
129       SMESH_subMesh* faceSM = mesh->GetMesh()->GetSubMesh( face );
130       faceSM->SetEventListener( this, new _Data( mesh ), faceSM );
131     }
132     // Finds a proxy mesh of face
133     static SMESH_ProxyMesh::Ptr FindProxyMeshOfFace( const TopoDS_Shape& face,
134                                                      SMESH_Mesh&         mesh )
135     {
136       SMESH_ProxyMesh::Ptr proxy;
137       SMESH_subMesh* faceSM = mesh.GetSubMesh( face );
138       if ( EventListenerData* ld = faceSM->GetEventListenerData( Name() ))
139         proxy = static_cast< _Data* >( ld )->_mesh;
140       return proxy;
141     }
142     // Treat events
143     void ProcessEvent(const int          event,
144                       const int          eventType,
145                       SMESH_subMesh*     subMesh,
146                       EventListenerData* data,
147                       const SMESH_Hypothesis*  /*hyp*/)
148     {
149       if ( event == SMESH_subMesh::CLEAN && eventType == SMESH_subMesh::COMPUTE_EVENT)
150         ((_Data*) data)->_mesh.reset();
151     }
152   private:
153     // holder of a proxy mesh
154     struct _Data : public SMESH_subMeshEventListenerData
155     {
156       SMESH_ProxyMesh::Ptr _mesh;
157       _Data( SMESH_ProxyMesh::Ptr& mesh )
158         :SMESH_subMeshEventListenerData( /*isDeletable=*/true), _mesh( mesh )
159       {}
160     };
161     // Returns identifier string
162     static const char* Name() { return "VISCOUS_2D::_ProxyMeshHolder"; }
163   };
164   
165   struct _PolyLine;
166   //--------------------------------------------------------------------------------
167   /*!
168    * \brief Segment connecting inner ends of two _LayerEdge's.
169    */
170   struct _Segment
171   {
172     const gp_XY* _uv[2];       // poiter to _LayerEdge::_uvIn
173     int          _indexInLine; // position in _PolyLine
174
175     _Segment() {}
176     _Segment(const gp_XY& p1, const gp_XY& p2):_indexInLine(-1) { _uv[0] = &p1; _uv[1] = &p2; }
177     const gp_XY& p1() const { return *_uv[0]; }
178     const gp_XY& p2() const { return *_uv[1]; }
179   };
180   //--------------------------------------------------------------------------------
181   /*!
182    * \brief Tree of _Segment's used for a faster search of _Segment's.
183    */
184   struct _SegmentTree : public SMESH_Quadtree
185   {
186     typedef boost::shared_ptr< _SegmentTree > Ptr;
187
188     _SegmentTree( const vector< _Segment >& segments );
189     void GetSegmentsNear( const _Segment& seg, vector< const _Segment* >& found );
190     void GetSegmentsNear( const gp_Ax2d& ray, vector< const _Segment* >& found );
191   protected:
192     _SegmentTree() {}
193     _SegmentTree* newChild() const { return new _SegmentTree; }
194     void          buildChildrenData();
195     Bnd_B2d*      buildRootBox();
196   private:
197     static int    maxNbSegInLeaf() { return 5; }
198     struct _SegBox
199     {
200       const _Segment* _seg;
201       bool            _iMin[2];
202       void Set( const _Segment& seg )
203       {
204         _seg = &seg;
205         _iMin[0] = ( seg._uv[1]->X() < seg._uv[0]->X() );
206         _iMin[1] = ( seg._uv[1]->Y() < seg._uv[0]->Y() );
207       }
208       bool IsOut( const _Segment& seg ) const;
209       bool IsOut( const gp_Ax2d& ray ) const;
210     };
211     vector< _SegBox > _segments;
212   };
213   //--------------------------------------------------------------------------------
214   /*!
215    * \brief Edge normal to FACE boundary, connecting a point on EDGE (_uvOut)
216    * and a point of a layer internal boundary (_uvIn)
217    */
218   struct _LayerEdge
219   {
220     gp_XY         _uvOut;    // UV on the FACE boundary
221     gp_XY         _uvIn;     // UV inside the FACE
222     double        _length2D; // distance between _uvOut and _uvIn
223
224     bool          _isBlocked;// is more inflation possible or not
225
226     gp_XY         _normal2D; // to pcurve
227     double        _len2dTo3dRatio; // to pass 2D <--> 3D
228     gp_Ax2d       _ray;      // a ray starting at _uvOut
229
230     vector<gp_XY> _uvRefined; // divisions by layers
231
232     bool SetNewLength( const double length );
233   };
234   //--------------------------------------------------------------------------------
235   /*!
236    * \brief Poly line composed of _Segment's of one EDGE.
237    *        It's used to detect intersection of inflated layers by intersecting
238    *        _Segment's in 2D.
239    */
240   struct _PolyLine
241   {
242     StdMeshers_FaceSide* _wire;
243     int                  _edgeInd;     // index of my EDGE in _wire
244     bool                 _advancable;  // true if there is a viscous layer on my EDGE
245     bool                 _isStraight2D;// pcurve type
246     _PolyLine*           _leftLine;    // lines of neighbour EDGE's
247     _PolyLine*           _rightLine;
248     int                  _firstPntInd; // index in vector<UVPtStruct> of _wire
249     int                  _lastPntInd;
250
251     vector< _LayerEdge > _lEdges;      /* _lEdges[0] is usually is not treated
252                                           as it is equal to the last one of the _leftLine */
253     vector< _Segment >   _segments;    // segments connecting _uvIn's of _lEdges
254     _SegmentTree::Ptr    _segTree;
255
256     vector< _PolyLine* > _reachableLines;       // lines able to interfere with my layer
257
258     vector< const SMDS_MeshNode* > _leftNodes;  // nodes built from a left VERTEX
259     vector< const SMDS_MeshNode* > _rightNodes; // nodes built from a right VERTEX
260
261     typedef vector< _Segment >::iterator   TSegIterator;
262     typedef vector< _LayerEdge >::iterator TEdgeIterator;
263
264     TIDSortedElemSet     _newFaces; // faces generated from this line
265
266     bool IsCommonEdgeShared( const _PolyLine& other );
267     size_t FirstLEdge() const
268     {
269       return ( _leftLine->_advancable && _lEdges.size() > 2 ) ? 1 : 0;
270     }
271     bool IsAdjacent( const _Segment& seg, const _LayerEdge* LE=0 ) const
272     {
273       if ( LE && seg._indexInLine < _lEdges.size() &&
274            ( seg._uv[0] == & LE->_uvIn ||
275              seg._uv[1] == & LE->_uvIn ))
276         return true;
277       return ( & seg == &_leftLine->_segments.back() ||
278                & seg == &_rightLine->_segments[0] );
279     }
280   };
281   //--------------------------------------------------------------------------------
282   /*!
283    * \brief Intersector of _Segment's
284    */
285   struct _SegmentIntersection
286   {
287     gp_XY    _vec1, _vec2;     // Vec( _seg.p1(), _seg.p2() )
288     gp_XY    _vec21;           // Vec( _seg2.p1(), _seg1.p1() )
289     double   _D;               // _vec1.Crossed( _vec2 )
290     double   _param1, _param2; // intersection param on _seg1 and _seg2
291
292     bool Compute(const _Segment& seg1, const _Segment& seg2, bool seg2IsRay = false )
293     {
294       const double eps = 1e-10;
295       _vec1  = seg1.p2() - seg1.p1(); 
296       _vec2  = seg2.p2() - seg2.p1(); 
297       _vec21 = seg1.p1() - seg2.p1(); 
298       _D = _vec1.Crossed(_vec2);
299       if ( fabs(_D) < std::numeric_limits<double>::min())
300         return false;
301       _param1 = _vec2.Crossed(_vec21) / _D; 
302       if (_param1 < -eps || _param1 > 1 + eps )
303         return false;
304       _param2 = _vec1.Crossed(_vec21) / _D; 
305       if (_param2 < -eps || ( !seg2IsRay && _param2 > 1 + eps ))
306         return false;
307       return true;
308     }
309     bool Compute( const _Segment& seg1, const gp_Ax2d& ray )
310     {
311       gp_XY segEnd = ray.Location().XY() + ray.Direction().XY();
312       _Segment seg2( ray.Location().XY(), segEnd );
313       return Compute( seg1, seg2, true );
314     }
315     //gp_XY GetPoint() { return _seg1.p1() + _param1 * _vec1; }
316   };
317   //--------------------------------------------------------------------------------
318
319   typedef map< const SMDS_MeshNode*, _LayerEdge*, TIDCompare > TNode2Edge;
320   
321   //--------------------------------------------------------------------------------
322   /*!
323    * \brief Builder of viscous layers
324    */
325   class _ViscousBuilder2D
326   {
327   public:
328     _ViscousBuilder2D(SMESH_Mesh&                       theMesh,
329                       const TopoDS_Face&                theFace,
330                       const StdMeshers_ViscousLayers2D* theHyp);
331     SMESH_ComputeErrorPtr GetError() const { return _error; }
332     // does it's job
333     SMESH_ProxyMesh::Ptr  Compute();
334
335   private:
336
337     bool findEdgesWithLayers();
338     bool makePolyLines();
339     bool inflate();
340     bool fixCollisions();
341     bool refine();
342     bool shrink();
343     bool improve();
344     bool toShrinkForAdjacent( const TopoDS_Face& adjFace,
345                               const TopoDS_Edge& E,
346                               const TopoDS_Vertex& V);
347     void setLenRatio( _LayerEdge& LE, const gp_Pnt& pOut );
348     void setLayerEdgeData( _LayerEdge&           lEdge,
349                            const double          u,
350                            Handle(Geom2d_Curve)& pcurve,
351                            const bool            reverse);
352     void adjustCommonEdge( _PolyLine& LL, _PolyLine& LR );
353     void calcLayersHeight(const double    totalThick,
354                           vector<double>& heights);
355     bool removeMeshFaces(const TopoDS_Shape& face);
356
357     bool              error( const string& text );
358     SMESHDS_Mesh*     getMeshDS() { return _mesh->GetMeshDS(); }
359     _ProxyMeshOfFace* getProxyMesh();
360
361     // debug
362     //void makeGroupOfLE();
363
364   private:
365
366     // input data
367     SMESH_Mesh*                 _mesh;
368     TopoDS_Face                 _face;
369     const StdMeshers_ViscousLayers2D* _hyp;
370
371     // result data
372     SMESH_ProxyMesh::Ptr        _proxyMesh;
373     SMESH_ComputeErrorPtr       _error;
374
375     // working data
376     Handle(Geom_Surface)        _surface;
377     SMESH_MesherHelper          _helper;
378     TSideVector                 _faceSideVec; // wires (StdMeshers_FaceSide) of _face
379     vector<_PolyLine>           _polyLineVec; // fronts to advance
380
381     double                      _fPowN; // to compute thickness of layers
382     double                      _thickness; // required or possible layers thickness
383
384     // sub-shapes of _face 
385     set<TGeomID>                _ignoreShapeIds; // ids of EDGEs w/o layers
386     set<TGeomID>                _noShrinkVert;   // ids of VERTEXes that are extremities
387     // of EDGEs along which _LayerEdge can't be inflated because no viscous layers
388     // defined on neighbour FACEs sharing an EDGE. Nonetheless _LayerEdge's
389     // are inflated along such EDGEs but then such _LayerEdge's are turned into
390     // a node on VERTEX, i.e. all nodes on a _LayerEdge are melded into one node.
391     
392   };
393
394   //================================================================================
395   /*!
396    * \brief Returns StdMeshers_ViscousLayers2D for the FACE
397    */
398   const StdMeshers_ViscousLayers2D* findHyp(SMESH_Mesh&        theMesh,
399                                             const TopoDS_Face& theFace)
400   {
401     SMESH_HypoFilter hypFilter
402       ( SMESH_HypoFilter::HasName( StdMeshers_ViscousLayers2D::GetHypType() ));
403     const SMESH_Hypothesis * hyp =
404       theMesh.GetHypothesis( theFace, hypFilter, /*ancestors=*/true );
405     return dynamic_cast< const StdMeshers_ViscousLayers2D* > ( hyp );
406   }
407
408 } // namespace VISCOUS_2D
409
410 //================================================================================
411 // StdMeshers_ViscousLayers hypothesis
412 //
413 StdMeshers_ViscousLayers2D::StdMeshers_ViscousLayers2D(int hypId, int studyId, SMESH_Gen* gen)
414   :StdMeshers_ViscousLayers(hypId, studyId, gen)
415 {
416   _name = StdMeshers_ViscousLayers2D::GetHypType();
417   _param_algo_dim = -2; // auxiliary hyp used by 2D algos
418 }
419 // --------------------------------------------------------------------------------
420 bool StdMeshers_ViscousLayers2D::SetParametersByMesh(const SMESH_Mesh*   theMesh,
421                                                      const TopoDS_Shape& theShape)
422 {
423   // TODO ???
424   return false;
425 }
426 // --------------------------------------------------------------------------------
427 SMESH_ProxyMesh::Ptr
428 StdMeshers_ViscousLayers2D::Compute(SMESH_Mesh&        theMesh,
429                                     const TopoDS_Face& theFace)
430 {
431   SMESH_ProxyMesh::Ptr pm;
432
433   const StdMeshers_ViscousLayers2D* vlHyp = VISCOUS_2D::findHyp( theMesh, theFace );
434   if ( vlHyp )
435   {
436     VISCOUS_2D::_ViscousBuilder2D builder( theMesh, theFace, vlHyp );
437     pm = builder.Compute();
438     SMESH_ComputeErrorPtr error = builder.GetError();
439     if ( error && !error->IsOK() )
440       theMesh.GetSubMesh( theFace )->GetComputeError() = error;
441     else if ( !pm )
442       pm.reset( new SMESH_ProxyMesh( theMesh ));
443     if ( getenv("ONLY_VL2D"))
444       pm.reset();
445   }
446   else
447   {
448     pm.reset( new SMESH_ProxyMesh( theMesh ));
449   }
450   return pm;
451 }
452 // --------------------------------------------------------------------------------
453 void StdMeshers_ViscousLayers2D::RestoreListeners() const
454 {
455   StudyContextStruct* sc = _gen->GetStudyContext( _studyId );
456   std::map < int, SMESH_Mesh * >::iterator i_smesh = sc->mapMesh.begin();
457   for ( ; i_smesh != sc->mapMesh.end(); ++i_smesh )
458   {
459     SMESH_Mesh* smesh = i_smesh->second;
460     if ( !smesh ||
461          !smesh->HasShapeToMesh() ||
462          !smesh->GetMeshDS() ||
463          !smesh->GetMeshDS()->IsUsedHypothesis( this ))
464       continue;
465
466     // set event listeners to EDGE's of FACE where this hyp is used
467     TopoDS_Shape shape = i_smesh->second->GetShapeToMesh();
468     for ( TopExp_Explorer face( shape, TopAbs_FACE); face.More(); face.Next() )
469       if ( SMESH_Algo* algo = _gen->GetAlgo( *smesh, face.Current() ))
470       {
471         const std::list <const SMESHDS_Hypothesis *> & usedHyps =
472           algo->GetUsedHypothesis( *smesh, face.Current(), /*ignoreAuxiliary=*/false );
473         if ( std::find( usedHyps.begin(), usedHyps.end(), this ) != usedHyps.end() )
474           for ( TopExp_Explorer edge( face.Current(), TopAbs_EDGE); edge.More(); edge.Next() )
475             VISCOUS_3D::ToClearSubWithMain( smesh->GetSubMesh( edge.Current() ), face.Current() );
476       }
477   }
478 }
479 // END StdMeshers_ViscousLayers2D hypothesis
480 //================================================================================
481
482 using namespace VISCOUS_2D;
483
484 //================================================================================
485 /*!
486  * \brief Constructor of _ViscousBuilder2D
487  */
488 //================================================================================
489
490 _ViscousBuilder2D::_ViscousBuilder2D(SMESH_Mesh&                       theMesh,
491                                      const TopoDS_Face&                theFace,
492                                      const StdMeshers_ViscousLayers2D* theHyp):
493   _mesh( &theMesh ), _face( theFace ), _hyp( theHyp ), _helper( theMesh )
494 {
495   _helper.SetSubShape( _face );
496   _helper.SetElementsOnShape(true);
497
498   _surface = BRep_Tool::Surface( theFace );
499
500   if ( _hyp )
501     _fPowN = pow( _hyp->GetStretchFactor(), _hyp->GetNumberLayers() );
502 }
503
504 //================================================================================
505 /*!
506  * \brief Stores error description and returns false
507  */
508 //================================================================================
509
510 bool _ViscousBuilder2D::error(const string& text )
511 {
512   cout << "_ViscousBuilder2D::error " << text << endl;
513   _error->myName    = COMPERR_ALGO_FAILED;
514   _error->myComment = string("Viscous layers builder 2D: ") + text;
515   if ( SMESH_subMesh* sm = _mesh->GetSubMesh( _face ) )
516   {
517     SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
518     if ( smError && smError->myAlgo )
519       _error->myAlgo = smError->myAlgo;
520     smError = _error;
521   }
522   //makeGroupOfLE(); // debug
523
524   return false;
525 }
526
527 //================================================================================
528 /*!
529  * \brief Does its job
530  */
531 //================================================================================
532
533 SMESH_ProxyMesh::Ptr _ViscousBuilder2D::Compute()
534 {
535   _error       = SMESH_ComputeError::New(COMPERR_OK);
536   _faceSideVec = StdMeshers_FaceSide::GetFaceWires( _face, *_mesh, true, _error );
537   if ( !_error->IsOK() )
538     return _proxyMesh;
539
540   if ( !findEdgesWithLayers() ) // analysis of a shape
541     return _proxyMesh;
542
543   if ( ! makePolyLines() ) // creation of fronts
544     return _proxyMesh;
545     
546   if ( ! inflate() ) // advance fronts
547     return _proxyMesh;
548
549   if ( !shrink() ) // shrink segments on edges w/o layers
550     return _proxyMesh;
551
552   if ( ! refine() ) // make faces
553     return _proxyMesh;
554
555   //improve();
556
557   return _proxyMesh;
558 }
559
560 //================================================================================
561 /*!
562  * \brief Finds EDGE's to make viscous layers on.
563  */
564 //================================================================================
565
566 bool _ViscousBuilder2D::findEdgesWithLayers()
567 {
568   // collect all EDGEs to ignore defined by hyp
569   int nbMyEdgesIgnored = 0;
570   vector<TGeomID> ids = _hyp->GetBndShapesToIgnore();
571   for ( size_t i = 0; i < ids.size(); ++i )
572   {
573     const TopoDS_Shape& s = getMeshDS()->IndexToShape( ids[i] );
574     if ( !s.IsNull() && s.ShapeType() == TopAbs_EDGE ) {
575       _ignoreShapeIds.insert( ids[i] );
576       nbMyEdgesIgnored += ( _helper.IsSubShape( s, _face ));
577     }
578   }
579
580   // check all EDGEs of the _face
581   int totalNbEdges = 0;
582   for ( size_t iWire = 0; iWire < _faceSideVec.size(); ++iWire )
583   {
584     StdMeshers_FaceSidePtr wire = _faceSideVec[ iWire ];
585     totalNbEdges += wire->NbEdges();
586     for ( int iE = 0; iE < wire->NbEdges(); ++iE )
587       if ( _helper.NbAncestors( wire->Edge( iE ), *_mesh, TopAbs_FACE ) > 1 )
588       {
589         // ignore internal EDGEs (shared by several FACEs)
590         TGeomID edgeID = getMeshDS()->ShapeToIndex( wire->Edge( iE ));
591         _ignoreShapeIds.insert( edgeID );
592
593         // check if ends of an EDGE are to be added to _noShrinkVert
594         PShapeIteratorPtr faceIt = _helper.GetAncestors( wire->Edge( iE ), *_mesh, TopAbs_FACE );
595         while ( const TopoDS_Shape* neighbourFace = faceIt->next() )
596         {
597           if ( neighbourFace->IsSame( _face )) continue;
598           SMESH_Algo* algo = _mesh->GetGen()->GetAlgo( *_mesh, *neighbourFace );
599           if ( !algo ) continue;
600
601           const StdMeshers_ViscousLayers2D* viscHyp = 0;
602           const list <const SMESHDS_Hypothesis *> & allHyps =
603             algo->GetUsedHypothesis(*_mesh, *neighbourFace, /*noAuxiliary=*/false);
604           list< const SMESHDS_Hypothesis *>::const_iterator hyp = allHyps.begin();
605           for ( ; hyp != allHyps.end() && !viscHyp; ++hyp )
606             viscHyp = dynamic_cast<const StdMeshers_ViscousLayers2D*>( *hyp );
607
608           set<TGeomID> neighbourIgnoreEdges;
609           if (viscHyp) {
610             vector<TGeomID> ids = _hyp->GetBndShapesToIgnore();
611             neighbourIgnoreEdges.insert( ids.begin(), ids.end() );
612           }
613           for ( int iV = 0; iV < 2; ++iV )
614           {
615             TopoDS_Vertex vertex = iV ? wire->LastVertex(iE) : wire->FirstVertex(iE);
616             if ( !viscHyp )
617               _noShrinkVert.insert( getMeshDS()->ShapeToIndex( vertex ));
618             else
619             {
620               PShapeIteratorPtr edgeIt = _helper.GetAncestors( vertex, *_mesh, TopAbs_EDGE );
621               while ( const TopoDS_Shape* edge = edgeIt->next() )
622                 if ( !edge->IsSame( wire->Edge( iE )) &&
623                      neighbourIgnoreEdges.count( getMeshDS()->ShapeToIndex( *edge )))
624                   _noShrinkVert.insert( getMeshDS()->ShapeToIndex( vertex ));
625             }
626           }
627         }
628       }
629   }
630   return ( nbMyEdgesIgnored < totalNbEdges );
631 }
632
633 //================================================================================
634 /*!
635  * \brief Create the inner front of the viscous layers and prepare data for infation
636  */
637 //================================================================================
638
639 bool _ViscousBuilder2D::makePolyLines()
640 {
641   // Create _PolyLines and _LayerEdge's
642
643   // count total nb of EDGEs to allocate _polyLineVec
644   int nbEdges = 0;
645   for ( size_t iWire = 0; iWire < _faceSideVec.size(); ++iWire )
646     nbEdges += _faceSideVec[ iWire ]->NbEdges();
647   _polyLineVec.resize( nbEdges );
648
649   // Assign data to _PolyLine's
650   // ---------------------------
651
652   size_t iPoLine = 0;
653   for ( size_t iWire = 0; iWire < _faceSideVec.size(); ++iWire )
654   {
655     StdMeshers_FaceSidePtr      wire = _faceSideVec[ iWire ];
656     const vector<UVPtStruct>& points = wire->GetUVPtStruct();
657     if ( points.empty() && wire->NbPoints() > 0 )
658       return error("Invalid node parameters on some EDGE");
659     int iPnt = 0;
660     for ( int iE = 0; iE < wire->NbEdges(); ++iE )
661     {
662       _PolyLine& L  = _polyLineVec[ iPoLine++ ];
663       L._wire       = wire.get();
664       L._edgeInd    = iE;
665       L._advancable = !_ignoreShapeIds.count( wire->EdgeID( iE ));
666
667       int iRight    = iPoLine - (( iE+1 < wire->NbEdges() ) ? 0 : wire->NbEdges() );
668       L._rightLine  = &_polyLineVec[ iRight ];
669       _polyLineVec[ iRight ]._leftLine = &L;
670
671       L._firstPntInd = iPnt;
672       double lastNormPar = wire->LastParameter( iE ) - 1e-10;
673       while ( points[ iPnt ].normParam < lastNormPar )
674         ++iPnt;
675       L._lastPntInd = iPnt;
676       L._lEdges.resize( Max( 3, L._lastPntInd - L._firstPntInd + 1 )); // 3 edges minimum
677
678       // TODO: add more _LayerEdge's to strongly curved EDGEs
679       // in order not to miss collisions
680
681       Handle(Geom2d_Curve) pcurve = L._wire->Curve2d( L._edgeInd );
682       const bool reverse = ( L._wire->Edge( iE ).Orientation() == TopAbs_REVERSED );
683       for ( int i = L._firstPntInd; i <= L._lastPntInd; ++i )
684       {
685         _LayerEdge& lEdge = L._lEdges[ i - L._firstPntInd ];
686         const double u = ( i == L._firstPntInd ? wire->FirstU(iE) : points[ i ].param );
687         setLayerEdgeData( lEdge, u, pcurve, reverse );
688         setLenRatio( lEdge, SMESH_TNodeXYZ( points[ i ].node ) );
689       }
690       if ( L._lastPntInd - L._firstPntInd + 1 < 3 ) // add 3d _LayerEdge in the middle
691       {
692         L._lEdges[2] = L._lEdges[1];
693         const double u = 0.5 * ( wire->FirstU(iE) + wire->LastU(iE) );
694         setLayerEdgeData( L._lEdges[1], u, pcurve, reverse );
695         gp_Pnt p = 0.5 * ( SMESH_TNodeXYZ( points[ L._firstPntInd ].node ) +
696                            SMESH_TNodeXYZ( points[ L._lastPntInd ].node ));
697         setLenRatio( L._lEdges[1], p );
698       }
699     }
700   }
701
702   // Fill _PolyLine's with _segments
703   // --------------------------------
704
705   double maxLen2dTo3dRatio = 0;
706   for ( iPoLine = 0; iPoLine < _polyLineVec.size(); ++iPoLine )
707   {
708     _PolyLine& L = _polyLineVec[ iPoLine ];
709     L._segments.resize( L._lEdges.size() - 1 );
710     for ( size_t i = 1; i < L._lEdges.size(); ++i )
711     {
712       _Segment & S   = L._segments[i-1];
713       S._uv[0]       = & L._lEdges[i-1]._uvIn;
714       S._uv[1]       = & L._lEdges[i  ]._uvIn;
715       S._indexInLine = i-1;
716       if ( maxLen2dTo3dRatio < L._lEdges[i]._len2dTo3dRatio )
717         maxLen2dTo3dRatio = L._lEdges[i]._len2dTo3dRatio;
718     }
719     // // connect _PolyLine's with segments, the 1st _LayerEdge of every _PolyLine
720     // // becomes not connected to any segment
721     // if ( L._leftLine->_advancable )
722     //   L._segments[0]._uv[0] = & L._leftLine->_lEdges.back()._uvIn;
723
724     L._segTree.reset( new _SegmentTree( L._segments ));
725   }
726
727   // Evaluate max possible _thickness if required layers thickness seems too high
728   // ----------------------------------------------------------------------------
729
730   _thickness = _hyp->GetTotalThickness();
731   _SegmentTree::box_type faceBndBox2D;
732   for ( iPoLine = 0; iPoLine < _polyLineVec.size(); ++iPoLine )
733     faceBndBox2D.Add( *_polyLineVec[ iPoLine]._segTree->getBox() );
734   double boxTol = 1e-3 * sqrt( faceBndBox2D.SquareExtent() );
735   //
736   if ( _thickness * maxLen2dTo3dRatio > sqrt( faceBndBox2D.SquareExtent() ) / 10 )
737   {
738     vector< const _Segment* > foundSegs;
739     double maxPossibleThick = 0;
740     _SegmentIntersection intersection;
741     for ( size_t iL1 = 0; iL1 < _polyLineVec.size(); ++iL1 )
742     {
743       _PolyLine& L1 = _polyLineVec[ iL1 ];
744       _SegmentTree::box_type boxL1 = * L1._segTree->getBox();
745       boxL1.Enlarge( boxTol );
746       for ( size_t iL2 = iL1+1; iL2 < _polyLineVec.size(); ++iL2 )
747       {
748         _PolyLine& L2 = _polyLineVec[ iL2 ];
749         _SegmentTree::box_type boxL2 = * L2._segTree->getBox();
750         boxL2.Enlarge( boxTol );
751         if ( boxL1.IsOut( boxL2 ))
752           continue;
753         for ( size_t iLE = 1; iLE < L1._lEdges.size(); ++iLE )
754         {
755           foundSegs.clear();
756           L2._segTree->GetSegmentsNear( L1._lEdges[iLE]._ray, foundSegs );
757           for ( size_t i = 0; i < foundSegs.size(); ++i )
758             if ( intersection.Compute( *foundSegs[i], L1._lEdges[iLE]._ray ))
759             {
760               double  distToL2 = intersection._param2 / L1._lEdges[iLE]._len2dTo3dRatio;
761               double psblThick = distToL2 / ( 1 + L1._advancable + L2._advancable );
762               if ( maxPossibleThick < psblThick )
763                 maxPossibleThick = psblThick;
764             }
765         }
766       }
767     }
768     _thickness = Min( _hyp->GetTotalThickness(), maxPossibleThick );
769   }
770
771   // Adjust _LayerEdge's at _PolyLine's extremities
772   // -----------------------------------------------
773
774   for ( iPoLine = 0; iPoLine < _polyLineVec.size(); ++iPoLine )
775   {
776     _PolyLine& LL = _polyLineVec[ iPoLine ];
777     _PolyLine& LR = *LL._rightLine;
778     adjustCommonEdge( LL, LR );
779   }
780   // recreate _segments if some _LayerEdge's have been removed by adjustCommonEdge()
781   for ( iPoLine = 0; iPoLine < _polyLineVec.size(); ++iPoLine )
782   {
783     _PolyLine& L = _polyLineVec[ iPoLine ];
784     // if ( L._segments.size() ==  L._lEdges.size() - 1 )
785     //   continue;
786     L._segments.resize( L._lEdges.size() - 1 );
787     for ( size_t i = 1; i < L._lEdges.size(); ++i )
788     {
789       _Segment & S   = L._segments[i-1];
790       S._uv[0]       = & L._lEdges[i-1]._uvIn;
791       S._uv[1]       = & L._lEdges[i  ]._uvIn;
792       S._indexInLine = i-1;
793     }
794     L._segTree.reset( new _SegmentTree( L._segments ));
795   }
796   // connect _PolyLine's with segments, the 1st _LayerEdge of every _PolyLine
797   // becomes not connected to any segment
798   for ( iPoLine = 0; iPoLine < _polyLineVec.size(); ++iPoLine )
799   {
800     _PolyLine& L = _polyLineVec[ iPoLine ];
801     if ( L._leftLine->_advancable )
802       L._segments[0]._uv[0] = & L._leftLine->_lEdges.back()._uvIn;
803   }
804
805   // Fill _reachableLines.
806   // ----------------------
807
808   // compute bnd boxes taking into account the layers total thickness
809   vector< _SegmentTree::box_type > lineBoxes( _polyLineVec.size() );
810   for ( iPoLine = 0; iPoLine < _polyLineVec.size(); ++iPoLine )
811   {
812     lineBoxes[ iPoLine ] = *_polyLineVec[ iPoLine ]._segTree->getBox();
813     if ( _polyLineVec[ iPoLine ]._advancable )
814       lineBoxes[ iPoLine ].Enlarge( maxLen2dTo3dRatio * _thickness * 2 );
815   }
816   // _reachableLines
817   for ( iPoLine = 0; iPoLine < _polyLineVec.size(); ++iPoLine )
818   {
819     _PolyLine& L1 = _polyLineVec[ iPoLine ];
820     for ( size_t iL2 = 0; iL2 < _polyLineVec.size(); ++iL2 )
821     {
822       _PolyLine& L2 = _polyLineVec[ iL2 ];
823       if ( iPoLine == iL2 || lineBoxes[ iPoLine ].IsOut( lineBoxes[ iL2 ]))
824         continue;
825       if ( !L1._advancable && ( L1._leftLine == &L2 || L1._rightLine == &L2 ))
826         continue;
827       // check reachability by _LayerEdge's
828       int iDelta = 1; //Max( 1, L1._lEdges.size() / 100 );
829       for ( size_t iLE = 1; iLE < L1._lEdges.size(); iLE += iDelta )
830       {
831         _LayerEdge& LE = L1._lEdges[iLE];
832         if ( !lineBoxes[ iL2 ].IsOut ( LE._uvOut,
833                                        LE._uvOut + LE._normal2D *_thickness * LE._len2dTo3dRatio ))
834         {
835           L1._reachableLines.push_back( & L2 );
836           break;
837         }
838       }
839     }
840     // add self to _reachableLines
841     Geom2dAdaptor_Curve pcurve( L1._wire->Curve2d( L1._edgeInd ));
842     L1._isStraight2D = ( pcurve.GetType() == GeomAbs_Line );
843     if ( !L1._isStraight2D )
844     {
845       // TODO: check carefully
846       L1._reachableLines.push_back( & L1 );
847     }
848   }
849
850   return true;
851 }
852
853 //================================================================================
854 /*!
855  * \brief adjust common _LayerEdge of two adjacent _PolyLine's
856  *  \param LL - left _PolyLine
857  *  \param LR - right _PolyLine
858  */
859 //================================================================================
860
861 void _ViscousBuilder2D::adjustCommonEdge( _PolyLine& LL, _PolyLine& LR )
862 {
863   int nbAdvancableL = LL._advancable + LR._advancable;
864   if ( nbAdvancableL == 0 )
865     return;
866
867   _LayerEdge& EL = LL._lEdges.back();
868   _LayerEdge& ER = LR._lEdges.front();
869   gp_XY normL    = EL._normal2D;
870   gp_XY normR    = ER._normal2D;
871   gp_XY tangL ( normL.Y(), -normL.X() );
872
873   // set common direction to a VERTEX _LayerEdge shared by two _PolyLine's
874   gp_XY normCommon = ( normL * int( LL._advancable ) +
875                        normR * int( LR._advancable )).Normalized();
876   EL._normal2D = normCommon;
877   EL._ray.SetLocation ( EL._uvOut );
878   EL._ray.SetDirection( EL._normal2D );
879   if ( nbAdvancableL == 1 ) { // _normal2D is true normal (not average)
880     EL._isBlocked = true; // prevent intersecting with _Segments of _advancable line
881     EL._length2D  = 0;
882   }
883   // update _LayerEdge::_len2dTo3dRatio according to a new direction
884   const vector<UVPtStruct>& points = LL._wire->GetUVPtStruct();
885   setLenRatio( EL, SMESH_TNodeXYZ( points[ LL._lastPntInd ].node ));
886
887   ER = EL;
888
889   const double dotNormTang = normR * tangL;
890   const bool    largeAngle = Abs( dotNormTang ) > 0.2;
891   if ( largeAngle ) // not 180 degrees
892   {
893     // recompute _len2dTo3dRatio to take into account angle between EDGEs
894     gp_Vec2d oldNorm( LL._advancable ? normL : normR );
895     double angleFactor  = 1. / Max( 0.3, Cos( oldNorm.Angle( normCommon )));
896     EL._len2dTo3dRatio *= angleFactor;
897     ER._len2dTo3dRatio  = EL._len2dTo3dRatio;
898
899     gp_XY normAvg = ( normL + normR ).Normalized(); // average normal at VERTEX
900
901     if ( dotNormTang < 0. ) // ---------------------------- CONVEX ANGLE
902     {
903       // Remove _LayerEdge's intersecting the normAvg to avoid collisions
904       // during inflate().
905       //
906       // find max length of the VERTEX based _LayerEdge whose direction is normAvg
907       double maxLen2D       = _thickness * EL._len2dTo3dRatio;
908       const gp_XY& pCommOut = ER._uvOut;
909       gp_XY        pCommIn  = pCommOut + normAvg * maxLen2D;
910       _Segment segCommon( pCommOut, pCommIn );
911       _SegmentIntersection intersection;
912       vector< const _Segment* > foundSegs;
913       for ( size_t iL1 = 0; iL1 < _polyLineVec.size(); ++iL1 )
914       {
915         _PolyLine& L1 = _polyLineVec[ iL1 ];
916         const _SegmentTree::box_type* boxL1 = L1._segTree->getBox();
917         if ( boxL1->IsOut ( pCommOut, pCommIn ))
918           continue;
919         for ( size_t iLE = 1; iLE < L1._lEdges.size(); ++iLE )
920         {
921           foundSegs.clear();
922           L1._segTree->GetSegmentsNear( segCommon, foundSegs );
923           for ( size_t i = 0; i < foundSegs.size(); ++i )
924             if ( intersection.Compute( *foundSegs[i], segCommon ) &&
925                  intersection._param2 > 1e-10 )
926             {
927               double len2D = intersection._param2 * maxLen2D / ( 2 + L1._advancable );
928               if ( len2D < maxLen2D ) {
929                 maxLen2D = len2D;
930                 pCommIn  = pCommOut + normAvg * maxLen2D; // here length of segCommon changes
931               }
932             }
933         }
934       }
935
936       // remove _LayerEdge's intersecting segCommon
937       for ( int isR = 0; isR < 2; ++isR ) // loop on [ LL, LR ]
938       {
939         _PolyLine&                 L = isR ? LR : LL;
940         _PolyLine::TEdgeIterator eIt = isR ? L._lEdges.begin()+1 : L._lEdges.end()-2;
941         int                      dIt = isR ? +1 : -1;
942         if ( nbAdvancableL == 1 && L._advancable && normL * normR > -0.01 )
943           continue;  // obtuse internal angle
944         // at least 3 _LayerEdge's should remain in a _PolyLine
945         if ( L._lEdges.size() < 4 ) continue;
946         size_t iLE = 1;
947         _SegmentIntersection lastIntersection;
948         for ( ; iLE < L._lEdges.size(); ++iLE, eIt += dIt )
949         {
950           gp_XY uvIn = eIt->_uvOut + eIt->_normal2D * _thickness * eIt->_len2dTo3dRatio;
951           _Segment segOfEdge( eIt->_uvOut, uvIn );
952           if ( !intersection.Compute( segCommon, segOfEdge ))
953             break;
954           lastIntersection._param1 = intersection._param1;
955           lastIntersection._param2 = intersection._param2;
956         }
957         if ( iLE >= L._lEdges.size () - 1 )
958         {
959           // all _LayerEdge's intersect the segCommon, limit inflation
960           // of remaining 2 _LayerEdge's
961           vector< _LayerEdge > newEdgeVec( Min( 3, L._lEdges.size() ));
962           newEdgeVec.front() = L._lEdges.front();
963           newEdgeVec.back()  = L._lEdges.back();
964           if ( newEdgeVec.size() == 3 )
965             newEdgeVec[1] = L._lEdges[ L._lEdges.size() / 2 ];
966           L._lEdges.swap( newEdgeVec );
967           if ( !isR ) std::swap( lastIntersection._param1 , lastIntersection._param2 );
968           L._lEdges.front()._len2dTo3dRatio *= lastIntersection._param1; // ??
969           L._lEdges.back ()._len2dTo3dRatio *= lastIntersection._param2;
970         }
971         else if ( iLE != 1 )
972         {
973           // eIt points to the _LayerEdge not intersecting with segCommon
974           if ( isR )
975             LR._lEdges.erase( LR._lEdges.begin()+1, eIt );
976           else
977             LL._lEdges.erase( eIt, --LL._lEdges.end() );
978         }
979       }
980     }
981     else // ------------------------------------------ CONCAVE ANGLE
982     {
983       if ( nbAdvancableL == 1 )
984       {
985         // make that the _LayerEdge at VERTEX is not shared by LL and LR:
986         // different normals is a sign that they are not shared
987         _LayerEdge& notSharedEdge = LL._advancable ? LR._lEdges[0] : LL._lEdges.back();
988         _LayerEdge&    sharedEdge = LR._advancable ? LR._lEdges[0] : LL._lEdges.back();
989
990         notSharedEdge._normal2D.SetCoord( 0.,0. );
991         sharedEdge._normal2D     = normAvg;
992         sharedEdge._isBlocked    = false;
993         notSharedEdge._isBlocked = true;
994       }
995     }
996   }
997 }
998
999 //================================================================================
1000 /*!
1001  * \brief initialize data of a _LayerEdge
1002  */
1003 //================================================================================
1004
1005 void _ViscousBuilder2D::setLayerEdgeData( _LayerEdge&           lEdge,
1006                                           const double          u,
1007                                           Handle(Geom2d_Curve)& pcurve,
1008                                           const bool            reverse)
1009 {
1010   gp_Pnt2d uv; gp_Vec2d tangent;
1011   pcurve->D1( u, uv, tangent );
1012   tangent.Normalize();
1013   if ( reverse )
1014     tangent.Reverse();
1015   lEdge._uvOut = lEdge._uvIn = uv.XY();
1016   lEdge._normal2D.SetCoord( -tangent.Y(), tangent.X() );
1017   lEdge._ray.SetLocation( lEdge._uvOut );
1018   lEdge._ray.SetDirection( lEdge._normal2D );
1019   lEdge._isBlocked = false;
1020   lEdge._length2D  = 0;
1021 }
1022
1023 //================================================================================
1024 /*!
1025  * \brief Compute and set _LayerEdge::_len2dTo3dRatio
1026  */
1027 //================================================================================
1028
1029 void _ViscousBuilder2D::setLenRatio( _LayerEdge& LE, const gp_Pnt& pOut )
1030 {
1031   const double probeLen2d = 1e-3;
1032
1033   gp_Pnt2d p2d = LE._uvOut + LE._normal2D * probeLen2d;
1034   gp_Pnt   p3d = _surface->Value( p2d.X(), p2d.Y() );
1035   double len3d = p3d.Distance( pOut );
1036   if ( len3d < std::numeric_limits<double>::min() )
1037     LE._len2dTo3dRatio = std::numeric_limits<double>::min();
1038   else
1039     LE._len2dTo3dRatio = probeLen2d / len3d;
1040 }
1041
1042 //================================================================================
1043 /*!
1044  * \brief Increase length of _LayerEdge's to reach the required thickness of layers
1045  */
1046 //================================================================================
1047
1048 bool _ViscousBuilder2D::inflate()
1049 {
1050   // Limit size of inflation step by geometry size found by
1051   // itersecting _LayerEdge's with _Segment's
1052   double minSize = _thickness, maxSize = 0;
1053   vector< const _Segment* > foundSegs;
1054   _SegmentIntersection intersection;
1055   for ( size_t iL1 = 0; iL1 < _polyLineVec.size(); ++iL1 )
1056   {
1057     _PolyLine& L1 = _polyLineVec[ iL1 ];
1058     for ( size_t iL2 = 0; iL2 < L1._reachableLines.size(); ++iL2 )
1059     {
1060       _PolyLine& L2 = * L1._reachableLines[ iL2 ];
1061       for ( size_t iLE = 1; iLE < L1._lEdges.size(); ++iLE )
1062       {
1063         foundSegs.clear();
1064         L2._segTree->GetSegmentsNear( L1._lEdges[iLE]._ray, foundSegs );
1065         for ( size_t i = 0; i < foundSegs.size(); ++i )
1066           if ( ! L1.IsAdjacent( *foundSegs[i], & L1._lEdges[iLE] ) &&
1067                intersection.Compute( *foundSegs[i], L1._lEdges[iLE]._ray ))
1068           {
1069             double distToL2 = intersection._param2 / L1._lEdges[iLE]._len2dTo3dRatio;
1070             double     size = distToL2 / ( 1 + L1._advancable + L2._advancable );
1071             if ( size < minSize )
1072               minSize = size;
1073             if ( size > maxSize )
1074               maxSize = size;
1075           }
1076       }
1077     }
1078   }
1079   if ( minSize > maxSize ) // no collisions possible
1080     maxSize = _thickness;
1081 #ifdef __myDEBUG
1082   cout << "-- minSize = " << minSize << ", maxSize = " << maxSize << endl;
1083 #endif
1084
1085   double curThick = 0, stepSize = minSize;
1086   int nbSteps = 0;
1087   if ( maxSize > _thickness )
1088     maxSize = _thickness;
1089   while ( curThick < maxSize )
1090   {
1091     curThick += stepSize * 1.25;
1092     if ( curThick > _thickness )
1093       curThick = _thickness;
1094
1095     // Elongate _LayerEdge's
1096     for ( size_t iL = 0; iL < _polyLineVec.size(); ++iL )
1097     {
1098       _PolyLine& L = _polyLineVec[ iL ];
1099       if ( !L._advancable ) continue;
1100       bool lenChange = false;
1101       for ( size_t iLE = L.FirstLEdge(); iLE < L._lEdges.size(); ++iLE )
1102         lenChange |= L._lEdges[iLE].SetNewLength( curThick );
1103       // for ( int k=0; k<L._segments.size(); ++k)
1104       //   cout << "( " << L._segments[k].p1().X() << ", " <<L._segments[k].p1().Y() << " ) "
1105       //        << "( " << L._segments[k].p2().X() << ", " <<L._segments[k].p2().Y() << " ) "
1106       //        << endl;
1107       if ( lenChange )
1108         L._segTree.reset( new _SegmentTree( L._segments ));
1109     }
1110
1111     // Avoid intersection of _Segment's
1112     bool allBlocked = fixCollisions();
1113     if ( allBlocked )
1114     {
1115       break; // no more inflating possible
1116     }
1117     stepSize = Max( stepSize , _thickness / 10. );
1118     nbSteps++;
1119   }
1120
1121   // if (nbSteps == 0 )
1122   //   return error("failed at the very first inflation step");
1123
1124
1125   // remove _LayerEdge's of one line intersecting with each other
1126   for ( size_t iL = 0; iL < _polyLineVec.size(); ++iL )
1127   {
1128     _PolyLine& L = _polyLineVec[ iL ];
1129     if ( !L._advancable ) continue;
1130
1131     // replace an inactive (1st) _LayerEdge with an active one of a neighbour _PolyLine
1132     if ( /*!L._leftLine->_advancable &&*/ L.IsCommonEdgeShared( *L._leftLine ) ) {
1133       L._lEdges[0] = L._leftLine->_lEdges.back();
1134     }
1135     if ( !L._rightLine->_advancable && L.IsCommonEdgeShared( *L._rightLine ) ) {
1136       L._lEdges.back() = L._rightLine->_lEdges[0];
1137     }
1138
1139     _SegmentIntersection intersection;
1140     for ( int isR = 0; ( isR < 2 && L._lEdges.size() > 2 ); ++isR )
1141     {
1142       int nbRemove = 0, deltaIt = isR ? -1 : +1;
1143       _PolyLine::TEdgeIterator eIt = isR ? L._lEdges.end()-1 : L._lEdges.begin();
1144       if ( eIt->_length2D == 0 ) continue;
1145       _Segment seg1( eIt->_uvOut, eIt->_uvIn );
1146       for ( eIt += deltaIt; nbRemove < L._lEdges.size()-1; eIt += deltaIt )
1147       {
1148         _Segment seg2( eIt->_uvOut, eIt->_uvIn );
1149         if ( !intersection.Compute( seg1, seg2 ))
1150           break;
1151         ++nbRemove;
1152       }
1153       if ( nbRemove > 0 ) {
1154         if ( nbRemove == L._lEdges.size()-1 ) // 1st and last _LayerEdge's intersect
1155         {
1156           --nbRemove;
1157           _LayerEdge& L0 = L._lEdges.front();
1158           _LayerEdge& L1 = L._lEdges.back();
1159           L0._length2D *= intersection._param1 * 0.5;
1160           L1._length2D *= intersection._param2 * 0.5;
1161           L0._uvIn = L0._uvOut + L0._normal2D * L0._length2D;
1162           L1._uvIn = L1._uvOut + L1._normal2D * L1._length2D;
1163           if ( L.IsCommonEdgeShared( *L._leftLine ))
1164             L._leftLine->_lEdges.back() = L0;
1165         }
1166         if ( isR )
1167           L._lEdges.erase( L._lEdges.end()-nbRemove-1,
1168                            L._lEdges.end()-nbRemove );
1169         else
1170           L._lEdges.erase( L._lEdges.begin()+1,
1171                            L._lEdges.begin()+1+nbRemove );
1172       }
1173     }
1174   }
1175   return true;
1176 }
1177
1178 //================================================================================
1179 /*!
1180  * \brief Remove intersection of _PolyLine's
1181  */
1182 //================================================================================
1183
1184 bool _ViscousBuilder2D::fixCollisions()
1185 {
1186   // look for intersections of _Segment's by intersecting _LayerEdge's with
1187   // _Segment's
1188   //double maxStep = 0, minStep = 1e+100;
1189   vector< const _Segment* > foundSegs;
1190   _SegmentIntersection intersection;
1191
1192   list< pair< _LayerEdge*, double > > edgeLenLimitList;
1193   list< _LayerEdge* >                 blockedEdgesList;
1194
1195   for ( size_t iL1 = 0; iL1 < _polyLineVec.size(); ++iL1 )
1196   {
1197     _PolyLine& L1 = _polyLineVec[ iL1 ];
1198     //if ( !L1._advancable ) continue;
1199     for ( size_t iL2 = 0; iL2 < L1._reachableLines.size(); ++iL2 )
1200     {
1201       _PolyLine& L2 = * L1._reachableLines[ iL2 ];
1202       for ( size_t iLE = L1.FirstLEdge(); iLE < L1._lEdges.size(); ++iLE )
1203       {
1204         _LayerEdge& LE1 = L1._lEdges[iLE];
1205         if ( LE1._isBlocked ) continue;
1206         foundSegs.clear();
1207         L2._segTree->GetSegmentsNear( LE1._ray, foundSegs );
1208         for ( size_t i = 0; i < foundSegs.size(); ++i )
1209         {
1210           if ( ! L1.IsAdjacent( *foundSegs[i], &LE1 ) &&
1211                intersection.Compute( *foundSegs[i], LE1._ray ))
1212           {
1213             const double dist2DToL2 = intersection._param2;
1214             double         newLen2D = dist2DToL2 / 2;
1215             if ( newLen2D < 1.1 * LE1._length2D ) // collision!
1216             {
1217               if ( newLen2D < LE1._length2D )
1218               {
1219                 blockedEdgesList.push_back( &LE1 );
1220                 if ( L1._advancable )
1221                 {
1222                   edgeLenLimitList.push_back( make_pair( &LE1, newLen2D ));
1223                   blockedEdgesList.push_back( &L2._lEdges[ foundSegs[i]->_indexInLine     ]);
1224                   blockedEdgesList.push_back( &L2._lEdges[ foundSegs[i]->_indexInLine + 1 ]);
1225                 }
1226                 else // here dist2DToL2 < 0 and LE1._length2D == 0
1227                 {
1228                   _LayerEdge LE2[2] = { L2._lEdges[ foundSegs[i]->_indexInLine     ],
1229                                         L2._lEdges[ foundSegs[i]->_indexInLine + 1 ] };
1230                   _Segment outSeg2( LE2[0]._uvOut, LE2[1]._uvOut );
1231                   intersection.Compute( outSeg2, LE1._ray );
1232                   newLen2D = intersection._param2 / 2;
1233
1234                   edgeLenLimitList.push_back( make_pair( &LE2[0], newLen2D ));
1235                   edgeLenLimitList.push_back( make_pair( &LE2[1], newLen2D ));
1236                 }
1237               }
1238             }
1239           }
1240         }
1241       }
1242     }
1243   }
1244
1245   // set limited length to _LayerEdge's
1246   list< pair< _LayerEdge*, double > >::iterator edge2Len = edgeLenLimitList.begin();
1247   for ( ; edge2Len != edgeLenLimitList.end(); ++edge2Len )
1248   {
1249     _LayerEdge* LE = edge2Len->first;
1250     LE->SetNewLength( edge2Len->second / LE->_len2dTo3dRatio );
1251     LE->_isBlocked = true;
1252   }
1253
1254   // block inflation of _LayerEdge's
1255   list< _LayerEdge* >::iterator edge = blockedEdgesList.begin();
1256   for ( ; edge != blockedEdgesList.end(); ++edge )
1257     (*edge)->_isBlocked = true;
1258
1259   // find a not blocked _LayerEdge
1260   for ( size_t iL = 0; iL < _polyLineVec.size(); ++iL )
1261   {
1262     _PolyLine& L = _polyLineVec[ iL ];
1263     if ( !L._advancable ) continue;
1264     for ( size_t iLE = L.FirstLEdge(); iLE < L._lEdges.size(); ++iLE )
1265       if ( !L._lEdges[ iLE ]._isBlocked )
1266         return false;
1267   }
1268
1269   return true;
1270 }
1271
1272 //================================================================================
1273 /*!
1274  * \brief Create new edges and shrink edges existing on a non-advancable _PolyLine
1275  *        adjacent to an advancable one.
1276  */
1277 //================================================================================
1278
1279 bool _ViscousBuilder2D::shrink()
1280 {
1281   gp_Pnt2d uv; //gp_Vec2d tangent;
1282   _SegmentIntersection intersection;
1283   double sign;
1284
1285   for ( size_t iL1 = 0; iL1 < _polyLineVec.size(); ++iL1 )
1286   {
1287     _PolyLine& L = _polyLineVec[ iL1 ]; // line with no layers
1288     if ( L._advancable )
1289       continue;
1290     const int nbAdvancable = ( L._rightLine->_advancable + L._leftLine->_advancable );
1291     if ( nbAdvancable == 0 )
1292       continue;
1293
1294     const TopoDS_Edge&        E = L._wire->Edge      ( L._edgeInd );
1295     const int            edgeID = L._wire->EdgeID    ( L._edgeInd );
1296     const double        edgeLen = L._wire->EdgeLength( L._edgeInd );
1297     Handle(Geom2d_Curve) pcurve = L._wire->Curve2d   ( L._edgeInd );
1298     const bool     edgeReversed = ( E.Orientation() == TopAbs_REVERSED );
1299
1300     SMESH_MesherHelper helper( *_mesh ); // to create nodes and edges on E
1301     helper.SetSubShape( E );
1302     helper.SetElementsOnShape( true );
1303
1304     // Check a FACE adjacent to _face by E
1305     bool existingNodesFound = false;
1306     TopoDS_Face adjFace;
1307     PShapeIteratorPtr faceIt = _helper.GetAncestors( E, *_mesh, TopAbs_FACE );
1308     while ( const TopoDS_Shape* f = faceIt->next() )
1309       if ( !_face.IsSame( *f ))
1310       {
1311         adjFace = TopoDS::Face( *f );
1312         SMESH_ProxyMesh::Ptr pm = _ProxyMeshHolder::FindProxyMeshOfFace( adjFace, *_mesh );
1313         if ( !pm || pm->NbProxySubMeshes() == 0 )
1314         {
1315           // There are no viscous layers on an adjacent FACE, clear it's 2D mesh
1316           removeMeshFaces( adjFace );
1317         }
1318         else
1319         {
1320           // There are viscous layers on the adjacent FACE; shrink must be already done;
1321           //
1322           // copy layer nodes
1323           //
1324           const vector<UVPtStruct>& points = L._wire->GetUVPtStruct();
1325           int iPFrom = L._firstPntInd, iPTo = L._lastPntInd;
1326           if ( L._leftLine->_advancable )
1327           {
1328             vector<gp_XY>& uvVec = L._lEdges.front()._uvRefined;
1329             for ( int i = 0; i < _hyp->GetNumberLayers(); ++i ) {
1330               const UVPtStruct& uvPt = points[ iPFrom + i + 1 ];
1331               L._leftNodes.push_back( uvPt.node );
1332               uvVec.push_back ( pcurve->Value( uvPt.param ).XY() );
1333             }
1334           }
1335           if ( L._rightLine->_advancable )
1336           {
1337             vector<gp_XY>& uvVec = L._lEdges.back()._uvRefined;
1338             for ( int i = 0; i < _hyp->GetNumberLayers(); ++i ) {
1339               const UVPtStruct& uvPt = points[ iPTo - i - 1 ];
1340               L._rightNodes.push_back( uvPt.node );
1341               uvVec.push_back ( pcurve->Value( uvPt.param ).XY() );
1342             }
1343           }
1344           // make proxy sub-mesh data of present nodes
1345           //
1346           if ( L._leftLine->_advancable )  iPFrom += _hyp->GetNumberLayers();
1347           if ( L._rightLine->_advancable ) iPTo   -= _hyp->GetNumberLayers();
1348           UVPtStructVec nodeDataVec( & points[ iPFrom ], & points[ iPTo + 1 ]);
1349
1350           double normSize = nodeDataVec.back().normParam - nodeDataVec.front().normParam;
1351           for ( int iP = nodeDataVec.size()-1; iP >= 0 ; --iP )
1352             nodeDataVec[iP].normParam =
1353               ( nodeDataVec[iP].normParam - nodeDataVec[0].normParam ) / normSize;
1354
1355           const SMDS_MeshNode* n = nodeDataVec.front().node;
1356           if ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX )
1357             nodeDataVec.front().param = L._wire->FirstU( L._edgeInd );
1358           n = nodeDataVec.back().node;
1359           if ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX )
1360             nodeDataVec.back().param = L._wire->LastU( L._edgeInd );
1361
1362           _ProxyMeshOfFace::_EdgeSubMesh* myEdgeSM = getProxyMesh()->GetEdgeSubMesh( edgeID );
1363           myEdgeSM->SetUVPtStructVec( nodeDataVec );
1364
1365           existingNodesFound = true;
1366         }
1367       } // loop on FACEs sharing E
1368
1369     if ( existingNodesFound )
1370       continue; // nothing more to do in this case
1371
1372     double u1 = L._wire->FirstU( L._edgeInd ), uf = u1;
1373     double u2 = L._wire->LastU ( L._edgeInd ), ul = u2;
1374
1375     // a ratio to pass 2D <--> 1D
1376     const double len1D = 1e-3;
1377     const double len2D = pcurve->Value(uf).Distance( pcurve->Value(uf+len1D));
1378     double len1dTo2dRatio = len1D / len2D;
1379
1380     // Get length of existing segments (from an edge start to a node) and their nodes
1381     const vector<UVPtStruct>& points = L._wire->GetUVPtStruct();
1382     UVPtStructVec nodeDataVec( & points[ L._firstPntInd ],
1383                                & points[ L._lastPntInd + 1 ]);
1384     nodeDataVec.front().param = u1; // U on vertex is correct on only one of shared edges
1385     nodeDataVec.back ().param = u2;
1386     nodeDataVec.front().normParam = 0;
1387     nodeDataVec.back ().normParam = 1;
1388     vector< double > segLengths( nodeDataVec.size() - 1 );
1389     BRepAdaptor_Curve curve( E );
1390     for ( size_t iP = 1; iP < nodeDataVec.size(); ++iP )
1391     {
1392       const double len = GCPnts_AbscissaPoint::Length( curve, uf, nodeDataVec[iP].param );
1393       segLengths[ iP-1 ] = len;
1394     }
1395
1396     // Before
1397     //  n1    n2    n3    n4
1398     //  x-----x-----x-----x-----
1399     //  |  e1    e2    e3    e4
1400
1401     // After
1402     //  n1          n2    n3
1403     //  x-x-x-x-----x-----x----
1404     //  | | | |  e1    e2    e3
1405
1406     // Move first and last parameters on EDGE (U of n1) according to layers' thickness
1407     // and create nodes of layers on EDGE ( -x-x-x )
1408     int isRShrinkedForAdjacent;
1409     UVPtStructVec nodeDataForAdjacent;
1410     for ( int isR = 0; isR < 2; ++isR )
1411     {
1412       _PolyLine* L2 = isR ? L._rightLine : L._leftLine; // line with layers
1413       if ( !L2->_advancable &&
1414            !toShrinkForAdjacent( adjFace, E, L._wire->FirstVertex( L._edgeInd + isR )))
1415         continue;
1416
1417       double & u = isR ? u2 : u1; // param to move
1418       double  u0 = isR ? ul : uf; // init value of the param to move
1419       int  iPEnd = isR ? nodeDataVec.size() - 1 : 0;
1420
1421       _LayerEdge& nearLE = isR ? L._lEdges.back() : L._lEdges.front();
1422       _LayerEdge&  farLE = isR ? L._lEdges.front() : L._lEdges.back();
1423
1424       // try to find length of advancement along L by intersecting L with
1425       // an adjacent _Segment of L2
1426
1427       double& length2D = nearLE._length2D;
1428       double  length1D = 0;
1429       sign = ( isR ^ edgeReversed ) ? -1. : 1.;
1430       //pcurve->D1( u, uv, tangent );
1431
1432       bool isConvex = false;
1433       if ( L2->_advancable )
1434       {
1435         int iFSeg2 = isR ? 0 : L2->_segments.size() - 1;
1436         int iLSeg2 = isR ? 1 : L2->_segments.size() - 2;
1437         gp_XY uvLSeg2In  = L2->_lEdges[ iLSeg2 ]._uvIn;
1438         gp_XY uvLSeg2Out = L2->_lEdges[ iLSeg2 ]._uvOut;
1439         gp_XY uvFSeg2Out = L2->_lEdges[ iFSeg2 ]._uvOut;
1440         Handle(Geom2d_Line) seg2Line = new Geom2d_Line( uvLSeg2In, uvFSeg2Out - uvLSeg2Out );
1441
1442         Geom2dAdaptor_Curve edgeCurve( pcurve, Min( uf, ul ), Max( uf, ul ));
1443         Geom2dAdaptor_Curve seg2Curve( seg2Line );
1444         Geom2dInt_GInter     curveInt( edgeCurve, seg2Curve, 1e-7, 1e-7 );
1445         isConvex = ( curveInt.IsDone() && !curveInt.IsEmpty() );
1446         if ( isConvex ) {
1447           /*                   convex VERTEX */
1448           length1D = Abs( u - curveInt.Point( 1 ).ParamOnFirst() );
1449           double maxDist2d = 2 * L2->_lEdges[ iLSeg2 ]._length2D;
1450           isConvex = ( length1D < maxDist2d * len1dTo2dRatio );
1451                                                   /*  |L  seg2     
1452                                                    *  |  o---o--- 
1453                                                    *  | /    |    
1454                                                    *  |/     |  L2
1455                                                    *  x------x---      */
1456         }
1457         if ( !isConvex ) { /* concave VERTEX */   /*  o-----o--- 
1458                                                    *   \    |    
1459                                                    *    \   |  L2
1460                                                    *     x--x--- 
1461                                                    *    /        
1462                                                    * L /               */
1463           length2D = L2->_lEdges[ iFSeg2 ]._length2D;
1464           //if ( L2->_advancable ) continue;
1465         }
1466       }
1467       else // L2 is advancable but in the face adjacent by L
1468       {
1469         length2D = farLE._length2D;
1470         if ( length2D == 0 ) {
1471           _LayerEdge& neighborLE =
1472             ( isR ? L._leftLine->_lEdges.back() : L._rightLine->_lEdges.front() );
1473           length2D = neighborLE._length2D;
1474           if ( length2D == 0 )
1475             length2D = _thickness * nearLE._len2dTo3dRatio;
1476         }
1477       }
1478
1479       // move u to the internal boundary of layers
1480       double maxLen3D = Min( _thickness, edgeLen / ( 1 + nbAdvancable ));
1481       double maxLen2D = maxLen3D * nearLE._len2dTo3dRatio;
1482       if ( !length2D ) length2D = length1D / len1dTo2dRatio;
1483       if ( Abs( length2D ) > maxLen2D )
1484         length2D = maxLen2D;
1485       nearLE._uvIn = nearLE._uvOut + nearLE._normal2D * length2D;
1486
1487       u += length2D * len1dTo2dRatio * sign;
1488       nodeDataVec[ iPEnd ].param = u;
1489
1490       gp_Pnt2d newUV = pcurve->Value( u );
1491       nodeDataVec[ iPEnd ].u = newUV.X();
1492       nodeDataVec[ iPEnd ].v = newUV.Y();
1493
1494       // compute params of layers on L
1495       vector<double> heights;
1496       calcLayersHeight( u - u0, heights );
1497       //
1498       vector< double > params( heights.size() );
1499       for ( size_t i = 0; i < params.size(); ++i )
1500         params[ i ] = u0 + heights[ i ];
1501
1502       // create nodes of layers and edges between them
1503       vector< const SMDS_MeshNode* >& layersNode = isR ? L._rightNodes : L._leftNodes;
1504       vector<gp_XY>& nodeUV = ( isR ? L._lEdges.back() : L._lEdges[0] )._uvRefined;
1505       nodeUV.resize    ( _hyp->GetNumberLayers() );
1506       layersNode.resize( _hyp->GetNumberLayers() );
1507       const SMDS_MeshNode* vertexNode = nodeDataVec[ iPEnd ].node;
1508       const SMDS_MeshNode *  prevNode = vertexNode;
1509       for ( size_t i = 0; i < params.size(); ++i )
1510       {
1511         gp_Pnt p        = curve.Value( params[i] );
1512         layersNode[ i ] = helper.AddNode( p.X(), p.Y(), p.Z(), /*id=*/0, params[i] );
1513         nodeUV    [ i ] = pcurve->Value( params[i] ).XY();
1514         helper.AddEdge( prevNode, layersNode[ i ] );
1515         prevNode = layersNode[ i ];
1516       }
1517
1518       // store data of layer nodes made for adjacent FACE
1519       if ( !L2->_advancable )
1520       {
1521         isRShrinkedForAdjacent = isR;
1522         nodeDataForAdjacent.resize( _hyp->GetNumberLayers() );
1523
1524         size_t iFrw = 0, iRev = nodeDataForAdjacent.size()-1, *i = isR ? &iRev : &iFrw;
1525         nodeDataForAdjacent[ *i ] = points[ isR ? L._lastPntInd : L._firstPntInd ];
1526         nodeDataForAdjacent[ *i ].param     = u0;
1527         nodeDataForAdjacent[ *i ].normParam = isR;
1528         for ( ++iFrw, --iRev; iFrw < layersNode.size(); ++iFrw, --iRev )
1529         {
1530           nodeDataForAdjacent[ *i ].node  = layersNode[ iFrw - 1 ];
1531           nodeDataForAdjacent[ *i ].u     = nodeUV    [ iFrw - 1 ].X();
1532           nodeDataForAdjacent[ *i ].v     = nodeUV    [ iFrw - 1 ].Y();
1533           nodeDataForAdjacent[ *i ].param = params    [ iFrw - 1 ];
1534         }
1535       }
1536       // replace a node on vertex by a node of last (most internal) layer
1537       // in a segment on E
1538       SMDS_ElemIteratorPtr segIt = vertexNode->GetInverseElementIterator( SMDSAbs_Edge );
1539       const SMDS_MeshNode* segNodes[3];
1540       while ( segIt->more() )
1541       {
1542         const SMDS_MeshElement* segment = segIt->next();
1543         if ( segment->getshapeId() != edgeID ) continue;
1544         
1545         const int nbNodes = segment->NbNodes();
1546         for ( int i = 0; i < nbNodes; ++i )
1547         {
1548           const SMDS_MeshNode* n = segment->GetNode( i );
1549           segNodes[ i ] = ( n == vertexNode ? layersNode.back() : n );
1550         }
1551         getMeshDS()->ChangeElementNodes( segment, segNodes, nbNodes );
1552         break;
1553       }
1554       nodeDataVec[ iPEnd ].node = layersNode.back();
1555
1556     } // loop on the extremities of L
1557
1558     // Shrink edges to fit in between the layers at EDGE ends
1559
1560     const double newLength = GCPnts_AbscissaPoint::Length( curve, u1, u2 );
1561     const double lenRatio  = newLength / edgeLen * ( edgeReversed ? -1. : 1. );
1562     for ( size_t iP = 1; iP < nodeDataVec.size()-1; ++iP )
1563     {
1564       const SMDS_MeshNode* oldNode = nodeDataVec[iP].node;
1565
1566       GCPnts_AbscissaPoint discret( curve, segLengths[iP-1] * lenRatio, u1 );
1567       if ( !discret.IsDone() )
1568         throw SALOME_Exception(LOCALIZED("GCPnts_AbscissaPoint failed"));
1569
1570       nodeDataVec[iP].param = discret.Parameter();
1571       if ( oldNode->GetPosition()->GetTypeOfPosition() != SMDS_TOP_EDGE )
1572         throw SALOME_Exception(SMESH_Comment("ViscousBuilder2D: not SMDS_TOP_EDGE node position: ")
1573                                << oldNode->GetPosition()->GetTypeOfPosition()
1574                                << " of node " << oldNode->GetID());
1575       SMDS_EdgePosition* pos = static_cast<SMDS_EdgePosition*>( oldNode->GetPosition() );
1576       pos->SetUParameter( nodeDataVec[iP].param );
1577
1578       gp_Pnt newP = curve.Value( nodeDataVec[iP].param );
1579       getMeshDS()->MoveNode( oldNode, newP.X(), newP.Y(), newP.Z() );
1580
1581       gp_Pnt2d newUV = pcurve->Value( nodeDataVec[iP].param ).XY();
1582       nodeDataVec[iP].u         = newUV.X();
1583       nodeDataVec[iP].v         = newUV.Y();
1584       nodeDataVec[iP].normParam = segLengths[iP-1] / edgeLen;
1585       // nodeDataVec[iP].x         = segLengths[iP-1] / edgeLen;
1586       // nodeDataVec[iP].y         = segLengths[iP-1] / edgeLen;
1587     }
1588
1589     // add nodeDataForAdjacent to nodeDataVec
1590     if ( !nodeDataForAdjacent.empty() )
1591     {
1592       const double par1      = isRShrinkedForAdjacent ? u2 : uf;
1593       const double par2      = isRShrinkedForAdjacent ? ul : u1;
1594       const double shrinkLen = GCPnts_AbscissaPoint::Length( curve, par1, par2 );
1595
1596       // compute new normParam for nodeDataVec
1597       for ( size_t iP = 0; iP < nodeDataVec.size()-1; ++iP )
1598         nodeDataVec[iP+1].normParam = segLengths[iP] / ( edgeLen + shrinkLen );
1599       double normDelta = 1 - nodeDataVec.back().normParam;
1600       if ( !isRShrinkedForAdjacent )
1601         for ( size_t iP = 0; iP < nodeDataVec.size(); ++iP )
1602           nodeDataVec[iP].normParam += normDelta;
1603
1604       // compute new normParam for nodeDataForAdjacent
1605       const double deltaR = isRShrinkedForAdjacent ? nodeDataVec.back().normParam : 0;
1606       for ( size_t iP = !isRShrinkedForAdjacent; iP < nodeDataForAdjacent.size(); ++iP )
1607       {
1608         double lenFromPar1 =
1609           GCPnts_AbscissaPoint::Length( curve, par1, nodeDataForAdjacent[iP].param );
1610         nodeDataForAdjacent[iP].normParam = deltaR + normDelta * lenFromPar1 / shrinkLen;
1611       }
1612       // concatenate nodeDataVec and nodeDataForAdjacent
1613       nodeDataVec.insert( isRShrinkedForAdjacent ? nodeDataVec.end() : nodeDataVec.begin(),
1614                           nodeDataForAdjacent.begin(), nodeDataForAdjacent.end() );
1615     }
1616
1617     // create a proxy sub-mesh containing the moved nodes
1618     _ProxyMeshOfFace::_EdgeSubMesh* edgeSM = getProxyMesh()->GetEdgeSubMesh( edgeID );
1619     edgeSM->SetUVPtStructVec( nodeDataVec );
1620
1621     // set a sub-mesh event listener to remove just created edges when
1622     // "ViscousLayers2D" hypothesis is modified
1623     VISCOUS_3D::ToClearSubWithMain( _mesh->GetSubMesh( E ), _face );
1624
1625   } // loop on _polyLineVec
1626
1627   return true;
1628 }
1629
1630 //================================================================================
1631 /*!
1632  * \brief Returns true if there will be a shrinked mesh on EDGE E of FACE adjFace
1633  *        near VERTEX V
1634  */
1635 //================================================================================
1636
1637 bool _ViscousBuilder2D::toShrinkForAdjacent( const TopoDS_Face&   adjFace,
1638                                              const TopoDS_Edge&   E,
1639                                              const TopoDS_Vertex& V)
1640 {
1641   if ( const StdMeshers_ViscousLayers2D* vlHyp = findHyp( *_mesh, adjFace ))
1642   {
1643     VISCOUS_2D::_ViscousBuilder2D builder( *_mesh, adjFace, vlHyp );
1644     builder.findEdgesWithLayers();
1645
1646     PShapeIteratorPtr edgeIt = _helper.GetAncestors( V, *_mesh, TopAbs_EDGE );
1647     while ( const TopoDS_Shape* edgeAtV = edgeIt->next() )
1648     {
1649       if ( !edgeAtV->IsSame( E ) &&
1650            _helper.IsSubShape( *edgeAtV, adjFace ) &&
1651            !builder._ignoreShapeIds.count( getMeshDS()->ShapeToIndex( *edgeAtV )))
1652       {
1653         return true;
1654       }
1655     }
1656   }
1657   return false;
1658 }
1659   
1660 //================================================================================
1661 /*!
1662  * \brief Make faces
1663  */
1664 //================================================================================
1665
1666 bool _ViscousBuilder2D::refine()
1667 {
1668   // remove elements and nodes from _face
1669   removeMeshFaces( _face );
1670
1671   // store a proxyMesh in a sub-mesh
1672   // make faces on each _PolyLine
1673   vector< double > layersHeight;
1674   double prevLen2D = -1;
1675   for ( size_t iL = 0; iL < _polyLineVec.size(); ++iL )
1676   {
1677     _PolyLine& L = _polyLineVec[ iL ];
1678     if ( !L._advancable ) continue;
1679
1680     // replace an inactive (1st) _LayerEdge with an active one of a neighbour _PolyLine
1681     size_t iLE = 0, nbLE = L._lEdges.size();
1682     const bool leftEdgeShared  = L.IsCommonEdgeShared( *L._leftLine );
1683     const bool rightEdgeShared = L.IsCommonEdgeShared( *L._rightLine );
1684     if ( /*!L._leftLine->_advancable &&*/ leftEdgeShared )
1685     {
1686       L._lEdges[0] = L._leftLine->_lEdges.back();
1687       iLE += int( !L._leftLine->_advancable );
1688     }
1689     if ( !L._rightLine->_advancable && rightEdgeShared )
1690     {
1691       L._lEdges.back() = L._rightLine->_lEdges[0];
1692       --nbLE;
1693     }
1694
1695     // limit length of neighbour _LayerEdge's to avoid sharp change of layers thickness
1696     vector< double > segLen( L._lEdges.size() );
1697     segLen[0] = 0.0;
1698     for ( size_t i = 1; i < segLen.size(); ++i )
1699     {
1700       // accumulate length of segments
1701       double sLen = (L._lEdges[i-1]._uvOut - L._lEdges[i]._uvOut ).Modulus();
1702       segLen[i] = segLen[i-1] + sLen;
1703     }
1704     for ( int isR = 0; isR < 2; ++isR )
1705     {
1706       size_t iF = 0, iL = L._lEdges.size()-1;
1707       size_t *i = isR ? &iL : &iF;
1708       //size_t iRef = *i;
1709       _LayerEdge* prevLE = & L._lEdges[ *i ];
1710       double weight = 0;
1711       for ( ++iF, --iL; iF < L._lEdges.size()-1; ++iF, --iL )
1712       {
1713         _LayerEdge& LE = L._lEdges[*i];
1714         if ( prevLE->_length2D > 0 ) {
1715           gp_XY tangent ( LE._normal2D.Y(), -LE._normal2D.X() );
1716           weight += Abs( tangent * ( prevLE->_uvIn - LE._uvIn )) / segLen.back();
1717           gp_XY prevTang = ( LE._uvOut - prevLE->_uvOut );
1718           gp_XY prevNorm    = gp_XY( -prevTang.Y(), prevTang.X() );
1719           double prevProj   = prevNorm * ( prevLE->_uvIn - prevLE->_uvOut );
1720           if ( prevProj > 0 ) {
1721             prevProj /= prevTang.Modulus();
1722             if ( LE._length2D < prevProj )
1723               weight += 0.75 * ( 1 - weight ); // length decrease is more preferable
1724             LE._length2D  = weight * LE._length2D + ( 1 - weight ) * prevProj;
1725             LE._uvIn = LE._uvOut + LE._normal2D * LE._length2D;
1726           }
1727         }
1728         prevLE = & LE;
1729       }
1730     }
1731
1732     // calculate intermediate UV on _LayerEdge's ( _LayerEdge::_uvRefined )
1733     for ( ; iLE < nbLE; ++iLE )
1734     {
1735       _LayerEdge& LE = L._lEdges[iLE];
1736       if ( fabs( LE._length2D - prevLen2D ) > LE._length2D / 100. )
1737       {
1738         calcLayersHeight( LE._length2D, layersHeight );
1739         prevLen2D = LE._length2D;
1740       }
1741       for ( size_t i = 0; i < layersHeight.size(); ++i )
1742         LE._uvRefined.push_back( LE._uvOut + LE._normal2D * layersHeight[i] );
1743     }
1744
1745     // nodes to create 1 layer of faces
1746     vector< const SMDS_MeshNode* > outerNodes( L._lastPntInd - L._firstPntInd + 1 );
1747     vector< const SMDS_MeshNode* > innerNodes( L._lastPntInd - L._firstPntInd + 1 );
1748
1749     // initialize outerNodes by node on the L._wire
1750     const vector<UVPtStruct>& points = L._wire->GetUVPtStruct();
1751     for ( int i = L._firstPntInd; i <= L._lastPntInd; ++i )
1752       outerNodes[ i-L._firstPntInd ] = points[i].node;
1753
1754     // compute normalized [0;1] node parameters of outerNodes
1755     vector< double > normPar( L._lastPntInd - L._firstPntInd + 1 );
1756     const double
1757       normF    = L._wire->FirstParameter( L._edgeInd ),
1758       normL    = L._wire->LastParameter ( L._edgeInd ),
1759       normDist = normL - normF;
1760     for ( int i = L._firstPntInd; i <= L._lastPntInd; ++i )
1761       normPar[ i - L._firstPntInd ] = ( points[i].normParam - normF ) / normDist;
1762
1763     // Create layers of faces
1764
1765     int hasLeftNode  = ( !L._leftLine->_rightNodes.empty() && leftEdgeShared );
1766     int hasRightNode = ( !L._rightLine->_leftNodes.empty() && rightEdgeShared );
1767     size_t iS, iN0 = hasLeftNode, nbN = innerNodes.size() - hasRightNode;
1768     L._leftNodes .resize( _hyp->GetNumberLayers() );
1769     L._rightNodes.resize( _hyp->GetNumberLayers() );
1770     for ( int iF = 0; iF < _hyp->GetNumberLayers(); ++iF ) // loop on layers of faces
1771     {
1772       // get accumulated length of intermediate segments
1773       for ( iS = 1; iS < segLen.size(); ++iS )
1774       {
1775         double sLen = (L._lEdges[iS-1]._uvRefined[iF] - L._lEdges[iS]._uvRefined[iF] ).Modulus();
1776         segLen[iS] = segLen[iS-1] + sLen;
1777       }
1778       // normalize the accumulated length
1779       for ( iS = 1; iS < segLen.size(); ++iS )
1780         segLen[iS] /= segLen.back();
1781
1782       // create innerNodes
1783       iS = 0;
1784       for ( size_t i = iN0; i < nbN; ++i )
1785       {
1786         while ( normPar[i] > segLen[iS+1] )
1787           ++iS;
1788         double r = ( normPar[i] - segLen[iS] ) / ( segLen[iS+1] - segLen[iS] );
1789         gp_XY uv = r * L._lEdges[iS+1]._uvRefined[iF] + (1-r) * L._lEdges[iS]._uvRefined[iF];
1790         gp_Pnt p = _surface->Value( uv.X(), uv.Y() );
1791         innerNodes[i] = _helper.AddNode( p.X(), p.Y(), p.Z(), /*id=*/0, uv.X(), uv.Y() );
1792       }
1793       if ( hasLeftNode ) innerNodes.front() = L._leftLine->_rightNodes[ iF ];
1794       if ( hasRightNode ) innerNodes.back() = L._rightLine->_leftNodes[ iF ];
1795       L._rightNodes[ iF ] = innerNodes.back();
1796       L._leftNodes [ iF ] = innerNodes.front();
1797
1798       // create faces
1799       // TODO care of orientation
1800       for ( size_t i = 1; i < innerNodes.size(); ++i )
1801         if ( SMDS_MeshElement* f = _helper.AddFace( outerNodes[ i-1 ], outerNodes[ i ],
1802                                                     innerNodes[ i ],   innerNodes[ i-1 ]))
1803           L._newFaces.insert( L._newFaces.end(), f );
1804
1805       outerNodes.swap( innerNodes );
1806     }
1807     // faces between not shared _LayerEdge's (at concave VERTEX)
1808     for ( int isR = 0; isR < 2; ++isR )
1809     {
1810       if ( isR ? rightEdgeShared : leftEdgeShared)
1811         continue;
1812       vector< const SMDS_MeshNode* > &
1813         lNodes = (isR ? L._rightNodes : L._leftLine->_rightNodes ),
1814         rNodes = (isR ? L._rightLine->_leftNodes : L._leftNodes );
1815       if ( lNodes.empty() || rNodes.empty() || lNodes.size() != rNodes.size() )
1816         continue;
1817
1818       for ( size_t i = 1; i < lNodes.size(); ++i )
1819         _helper.AddFace( lNodes[ i-1 ], rNodes[ i-1 ],
1820                          rNodes[ i ],   lNodes[ i ]);
1821
1822       const UVPtStruct& ptOnVertex = points[ isR ? L._lastPntInd : L._firstPntInd ];
1823       _helper.AddFace( ptOnVertex.node, rNodes[ 0 ], lNodes[ 0 ]);
1824
1825       // update nodeDataVec of an adjacent _PolyLine
1826       // int iAdjEdge = isR ? L._rightLine->_edgeInd : L._leftLine->_edgeInd;
1827       // _ProxyMeshOfFace::_EdgeSubMesh* adjEdgeSM
1828       //   = getProxyMesh()->GetEdgeSubMesh( L._wire->EdgeID( iAdjEdge ));
1829       // const UVPtStructVec& nodeDataVec = adjEdgeSM->GetUVPtStructVec();
1830       // if ( !nodeDataVec.empty() )
1831       // {
1832       //   UVPtStruct ptOnVertex;
1833       //   _LayerEdge& LE = isR ? L._lEdges.back() : L._lEdges.front();
1834       //   ptOnVertex.u         = LE._uvRefined.back().X();
1835       //   ptOnVertex.v         = LE._uvRefined.back().Y();
1836       //   ptOnVertex.node      = isR ? L._rightNodes.back() : L._leftNodes.back();
1837       //   ptOnVertex.param     = isR ? L._wire->FirstU( iAdjEdge ) :L._wire->LastU( iAdjEdge );
1838       //   ptOnVertex.normParam = isR ? 1 : 0;
1839       //   ptOnVertex.x         = ptOnVertex.normParam;
1840       //   ptOnVertex.y         = ptOnVertex.normParam;
1841
1842       //   int iN  = isR ? _hyp->GetNumberLayers() : 0;
1843       //   int nbN = nodeDataVec.size() - ( isR ? 0 : _hyp->GetNumberLayers() );
1844       //   UVPtStructVec newNodeData( nodeDataVec.begin() + iN,
1845       //                              nodeDataVec.begin() + nbN );
1846       //   newNodeData.insert( isR ? newNodeData.begin() : newNodeData.end(), ptOnVertex );
1847       //   adjEdgeSM->SetUVPtStructVec( newNodeData );
1848       // }
1849     }
1850
1851     // Fill the _ProxyMeshOfFace
1852
1853     UVPtStructVec nodeDataVec( outerNodes.size() ); // outerNodes swapped with innerNodes
1854     for ( size_t i = 0; i < outerNodes.size(); ++i )
1855     {
1856       gp_XY uv = _helper.GetNodeUV( _face, outerNodes[i] );
1857       nodeDataVec[i].u         = uv.X();
1858       nodeDataVec[i].v         = uv.Y();
1859       nodeDataVec[i].node      = outerNodes[i];
1860       nodeDataVec[i].param     = points [i + L._firstPntInd].param;
1861       nodeDataVec[i].normParam = normPar[i];
1862       nodeDataVec[i].x         = normPar[i];
1863       nodeDataVec[i].y         = normPar[i];
1864     }
1865     nodeDataVec.front().param = L._wire->FirstU( L._edgeInd );
1866     nodeDataVec.back() .param = L._wire->LastU ( L._edgeInd );
1867
1868     _ProxyMeshOfFace::_EdgeSubMesh* edgeSM
1869       = getProxyMesh()->GetEdgeSubMesh( L._wire->EdgeID( L._edgeInd ));
1870     edgeSM->SetUVPtStructVec( nodeDataVec );
1871
1872   } // loop on _PolyLine's
1873
1874   return true;
1875 }
1876
1877 //================================================================================
1878 /*!
1879  * \brief Improve quality of the created mesh elements
1880  */
1881 //================================================================================
1882
1883 bool _ViscousBuilder2D::improve()
1884 {
1885   if ( !_proxyMesh )
1886     return false;
1887
1888   // fixed nodes on EDGE's
1889   std::set<const SMDS_MeshNode*> fixedNodes;
1890   for ( size_t iWire = 0; iWire < _faceSideVec.size(); ++iWire )
1891   {
1892     StdMeshers_FaceSidePtr      wire = _faceSideVec[ iWire ];
1893     const vector<UVPtStruct>& points = wire->GetUVPtStruct();
1894     for ( size_t i = 0; i < points.size(); ++i )
1895       fixedNodes.insert( fixedNodes.end(), points[i].node );
1896   }
1897   // fixed proxy nodes
1898   for ( size_t iL = 0; iL < _polyLineVec.size(); ++iL )
1899   {
1900     _PolyLine&         L = _polyLineVec[ iL ];
1901     const TopoDS_Edge& E = L._wire->Edge( L._edgeInd );
1902     if ( const SMESH_ProxyMesh::SubMesh* sm = _proxyMesh->GetProxySubMesh( E ))
1903     {
1904       const UVPtStructVec& points = sm->GetUVPtStructVec();
1905       for ( size_t i = 0; i < points.size(); ++i )
1906         fixedNodes.insert( fixedNodes.end(), points[i].node );
1907     }
1908     for ( size_t i = 0; i < L._rightNodes.size(); ++i )
1909       fixedNodes.insert( fixedNodes.end(), L._rightNodes[i] );
1910   }
1911
1912   // smoothing
1913   SMESH_MeshEditor editor( _mesh );
1914   for ( size_t iL = 0; iL < _polyLineVec.size(); ++iL )
1915   {
1916     _PolyLine& L = _polyLineVec[ iL ];
1917     if ( L._isStraight2D ) continue;
1918     // SMESH_MeshEditor::SmoothMethod how =
1919     //   L._isStraight2D ? SMESH_MeshEditor::LAPLACIAN : SMESH_MeshEditor::CENTROIDAL;
1920     //editor.Smooth( L._newFaces, fixedNodes, how, /*nbIt = */3 );
1921     //editor.Smooth( L._newFaces, fixedNodes, SMESH_MeshEditor::LAPLACIAN, /*nbIt = */1 );
1922     editor.Smooth( L._newFaces, fixedNodes, SMESH_MeshEditor::CENTROIDAL, /*nbIt = */3 );
1923   }
1924   return true;
1925 }
1926
1927 //================================================================================
1928 /*!
1929  * \brief Remove elements and nodes from a face
1930  */
1931 //================================================================================
1932
1933 bool _ViscousBuilder2D::removeMeshFaces(const TopoDS_Shape& face)
1934 {
1935   // we don't use SMESH_subMesh::ComputeStateEngine() because of a listener
1936   // which clears EDGEs together with _face.
1937   bool thereWereElems = false;
1938   SMESH_subMesh* sm = _mesh->GetSubMesh( face );
1939   if ( SMESHDS_SubMesh* smDS = sm->GetSubMeshDS() )
1940   {
1941     SMDS_ElemIteratorPtr eIt = smDS->GetElements();
1942     thereWereElems = eIt->more();
1943     while ( eIt->more() ) getMeshDS()->RemoveFreeElement( eIt->next(), smDS );
1944     SMDS_NodeIteratorPtr nIt = smDS->GetNodes();
1945     while ( nIt->more() ) getMeshDS()->RemoveFreeNode( nIt->next(), smDS );
1946   }
1947   sm->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
1948
1949   return thereWereElems;
1950 }
1951
1952 //================================================================================
1953 /*!
1954  * \brief Creates a _ProxyMeshOfFace and store it in a sub-mesh of FACE
1955  */
1956 //================================================================================
1957
1958 _ProxyMeshOfFace* _ViscousBuilder2D::getProxyMesh()
1959 {
1960   if ( _proxyMesh.get() )
1961     return (_ProxyMeshOfFace*) _proxyMesh.get();
1962
1963   _ProxyMeshOfFace* proxyMeshOfFace = new _ProxyMeshOfFace( *_mesh );
1964   _proxyMesh.reset( proxyMeshOfFace );
1965   new _ProxyMeshHolder( _face, _proxyMesh );
1966
1967   return proxyMeshOfFace;
1968 }
1969
1970 //================================================================================
1971 /*!
1972  * \brief Calculate height of layers for the given thickness. Height is measured
1973  *        from the outer boundary
1974  */
1975 //================================================================================
1976
1977 void _ViscousBuilder2D::calcLayersHeight(const double    totalThick,
1978                                          vector<double>& heights)
1979 {
1980   heights.resize( _hyp->GetNumberLayers() );
1981   double h0;
1982   if ( _fPowN - 1 <= numeric_limits<double>::min() )
1983     h0 = totalThick / _hyp->GetNumberLayers();
1984   else
1985     h0 = totalThick * ( _hyp->GetStretchFactor() - 1 )/( _fPowN - 1 );
1986
1987   double hSum = 0, hi = h0;
1988   for ( int i = 0; i < _hyp->GetNumberLayers(); ++i )
1989   {
1990     hSum += hi;
1991     heights[ i ] = hSum;
1992     hi *= _hyp->GetStretchFactor();
1993   }
1994 }
1995
1996 //================================================================================
1997 /*!
1998  * \brief Elongate this _LayerEdge
1999  */
2000 //================================================================================
2001
2002 bool _LayerEdge::SetNewLength( const double length3D )
2003 {
2004   if ( _isBlocked ) return false;
2005
2006   //_uvInPrev = _uvIn;
2007   _length2D = length3D * _len2dTo3dRatio;
2008   _uvIn     = _uvOut + _normal2D * _length2D;
2009   return true;
2010 }
2011
2012 //================================================================================
2013 /*!
2014  * \brief Return true if _LayerEdge at a common VERTEX between EDGEs with
2015  *  and w/o layer is common to the both _PolyLine's. If this is true, nodes
2016  *  of this _LayerEdge are inflated along a _PolyLine w/o layer, else the nodes
2017  *  are inflated along _normal2D of _LayerEdge of EDGE with layer
2018  */
2019 //================================================================================
2020
2021 bool _PolyLine::IsCommonEdgeShared( const _PolyLine& other )
2022 {
2023   const double tol = 1e-30;
2024
2025   if ( & other == _leftLine )
2026     return _lEdges[0]._normal2D.IsEqual( _leftLine->_lEdges.back()._normal2D, tol );
2027
2028   if ( & other == _rightLine )
2029     return _lEdges.back()._normal2D.IsEqual( _rightLine->_lEdges[0]._normal2D, tol );
2030
2031   return false;
2032 }
2033
2034 //================================================================================
2035 /*!
2036  * \brief Constructor of SegmentTree
2037  */
2038 //================================================================================
2039
2040 _SegmentTree::_SegmentTree( const vector< _Segment >& segments ):
2041   SMESH_Quadtree()
2042 {
2043   _segments.resize( segments.size() );
2044   for ( size_t i = 0; i < segments.size(); ++i )
2045     _segments[i].Set( segments[i] );
2046
2047   compute();
2048 }
2049
2050 //================================================================================
2051 /*!
2052  * \brief Return the maximal bnd box
2053  */
2054 //================================================================================
2055
2056 _SegmentTree::box_type* _SegmentTree::buildRootBox()
2057 {
2058   _SegmentTree::box_type* box = new _SegmentTree::box_type;
2059   for ( size_t i = 0; i < _segments.size(); ++i )
2060   {
2061     box->Add( *_segments[i]._seg->_uv[0] );
2062     box->Add( *_segments[i]._seg->_uv[1] );
2063   }
2064   return box;
2065 }
2066
2067 //================================================================================
2068 /*!
2069  * \brief Redistrubute _segments among children
2070  */
2071 //================================================================================
2072
2073 void _SegmentTree::buildChildrenData()
2074 {
2075   for ( int i = 0; i < _segments.size(); ++i )
2076     for (int j = 0; j < nbChildren(); j++)
2077       if ( !myChildren[j]->getBox()->IsOut( *_segments[i]._seg->_uv[0],
2078                                             *_segments[i]._seg->_uv[1] ))
2079         ((_SegmentTree*)myChildren[j])->_segments.push_back( _segments[i]);
2080
2081   SMESHUtils::FreeVector( _segments ); // = _elements.clear() + free memory
2082
2083   for (int j = 0; j < nbChildren(); j++)
2084   {
2085     _SegmentTree* child = static_cast<_SegmentTree*>( myChildren[j]);
2086     child->myIsLeaf = ( child->_segments.size() <= maxNbSegInLeaf() );
2087   }
2088 }
2089
2090 //================================================================================
2091 /*!
2092  * \brief Return elements which can include the point
2093  */
2094 //================================================================================
2095
2096 void _SegmentTree::GetSegmentsNear( const _Segment&            seg,
2097                                     vector< const _Segment* >& found )
2098 {
2099   if ( getBox()->IsOut( *seg._uv[0], *seg._uv[1] ))
2100     return;
2101
2102   if ( isLeaf() )
2103   {
2104     for ( int i = 0; i < _segments.size(); ++i )
2105       if ( !_segments[i].IsOut( seg ))
2106         found.push_back( _segments[i]._seg );
2107   }
2108   else
2109   {
2110     for (int i = 0; i < nbChildren(); i++)
2111       ((_SegmentTree*) myChildren[i])->GetSegmentsNear( seg, found );
2112   }
2113 }
2114
2115
2116 //================================================================================
2117 /*!
2118  * \brief Return segments intersecting a ray
2119  */
2120 //================================================================================
2121
2122 void _SegmentTree::GetSegmentsNear( const gp_Ax2d&             ray,
2123                                     vector< const _Segment* >& found )
2124 {
2125   if ( getBox()->IsOut( ray ))
2126     return;
2127
2128   if ( isLeaf() )
2129   {
2130     for ( int i = 0; i < _segments.size(); ++i )
2131       if ( !_segments[i].IsOut( ray ))
2132         found.push_back( _segments[i]._seg );
2133   }
2134   else
2135   {
2136     for (int i = 0; i < nbChildren(); i++)
2137       ((_SegmentTree*) myChildren[i])->GetSegmentsNear( ray, found );
2138   }
2139 }
2140
2141 //================================================================================
2142 /*!
2143  * \brief Classify a _Segment
2144  */
2145 //================================================================================
2146
2147 bool _SegmentTree::_SegBox::IsOut( const _Segment& seg ) const
2148 {
2149   const double eps = std::numeric_limits<double>::min();
2150   for ( int iC = 0; iC < 2; ++iC )
2151   {
2152     if ( seg._uv[0]->Coord(iC+1) < _seg->_uv[ _iMin[iC]]->Coord(iC+1)+eps &&
2153          seg._uv[1]->Coord(iC+1) < _seg->_uv[ _iMin[iC]]->Coord(iC+1)+eps )
2154       return true;
2155     if ( seg._uv[0]->Coord(iC+1) > _seg->_uv[ 1-_iMin[iC]]->Coord(iC+1)-eps &&
2156          seg._uv[1]->Coord(iC+1) > _seg->_uv[ 1-_iMin[iC]]->Coord(iC+1)-eps )
2157       return true;
2158   }
2159   return false;
2160 }
2161
2162 //================================================================================
2163 /*!
2164  * \brief Classify a ray
2165  */
2166 //================================================================================
2167
2168 bool _SegmentTree::_SegBox::IsOut( const gp_Ax2d& ray ) const
2169 {
2170   double distBoxCenter2Ray =
2171     ray.Direction().XY() ^ ( ray.Location().XY() - 0.5 * (*_seg->_uv[0] + *_seg->_uv[1]));
2172
2173   double boxSectionDiam =
2174     Abs( ray.Direction().X() ) * ( _seg->_uv[1-_iMin[1]]->Y() - _seg->_uv[_iMin[1]]->Y() ) +
2175     Abs( ray.Direction().Y() ) * ( _seg->_uv[1-_iMin[0]]->X() - _seg->_uv[_iMin[0]]->X() );
2176
2177   return Abs( distBoxCenter2Ray ) > 0.5 * boxSectionDiam;
2178 }