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