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