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