Salome HOME
Update copyrights 2014.
[modules/smesh.git] / src / StdMeshers / StdMeshers_ViscousLayers.cxx
1 // Copyright (C) 2007-2014  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, or (at your option) any later version.
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_ViscousLayers.cxx
21 // Created   : Wed Dec  1 15:15:34 2010
22 // Author    : Edward AGAPOV (eap)
23
24 #include "StdMeshers_ViscousLayers.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_MeshAlgos.hxx"
41 #include "SMESH_MesherHelper.hxx"
42 #include "SMESH_ProxyMesh.hxx"
43 #include "SMESH_subMesh.hxx"
44 #include "SMESH_subMeshEventListener.hxx"
45 #include "StdMeshers_FaceSide.hxx"
46
47 #include <BRepAdaptor_Curve2d.hxx>
48 #include <BRep_Tool.hxx>
49 #include <Bnd_B2d.hxx>
50 #include <Bnd_B3d.hxx>
51 #include <ElCLib.hxx>
52 #include <GCPnts_AbscissaPoint.hxx>
53 #include <Geom2d_Circle.hxx>
54 #include <Geom2d_Line.hxx>
55 #include <Geom2d_TrimmedCurve.hxx>
56 #include <GeomAdaptor_Curve.hxx>
57 #include <Geom_Circle.hxx>
58 #include <Geom_Curve.hxx>
59 #include <Geom_Line.hxx>
60 #include <Geom_TrimmedCurve.hxx>
61 #include <Precision.hxx>
62 #include <Standard_ErrorHandler.hxx>
63 #include <TColStd_Array1OfReal.hxx>
64 #include <TopExp.hxx>
65 #include <TopExp_Explorer.hxx>
66 #include <TopTools_IndexedMapOfShape.hxx>
67 #include <TopTools_ListOfShape.hxx>
68 #include <TopTools_MapOfShape.hxx>
69 #include <TopoDS.hxx>
70 #include <TopoDS_Edge.hxx>
71 #include <TopoDS_Face.hxx>
72 #include <TopoDS_Vertex.hxx>
73 #include <gp_Ax1.hxx>
74 #include <gp_Vec.hxx>
75 #include <gp_XY.hxx>
76
77 #include <list>
78 #include <string>
79 #include <cmath>
80 #include <limits>
81
82 //#define __myDEBUG
83
84 using namespace std;
85
86 //================================================================================
87 namespace VISCOUS_3D
88 {
89   typedef int TGeomID;
90
91   enum UIndex { U_TGT = 1, U_SRC, LEN_TGT };
92
93   /*!
94    * \brief SMESH_ProxyMesh computed by _ViscousBuilder for a SOLID.
95    * It is stored in a SMESH_subMesh of the SOLID as SMESH_subMeshEventListenerData
96    */
97   struct _MeshOfSolid : public SMESH_ProxyMesh,
98                         public SMESH_subMeshEventListenerData
99   {
100     bool _n2nMapComputed;
101
102     _MeshOfSolid( SMESH_Mesh* mesh)
103       :SMESH_subMeshEventListenerData( /*isDeletable=*/true),_n2nMapComputed(false)
104     {
105       SMESH_ProxyMesh::setMesh( *mesh );
106     }
107
108     // returns submesh for a geom face
109     SMESH_ProxyMesh::SubMesh* getFaceSubM(const TopoDS_Face& F, bool create=false)
110     {
111       TGeomID i = SMESH_ProxyMesh::shapeIndex(F);
112       return create ? SMESH_ProxyMesh::getProxySubMesh(i) : findProxySubMesh(i);
113     }
114     void setNode2Node(const SMDS_MeshNode*                 srcNode,
115                       const SMDS_MeshNode*                 proxyNode,
116                       const SMESH_ProxyMesh::SubMesh* subMesh)
117     {
118       SMESH_ProxyMesh::setNode2Node( srcNode,proxyNode,subMesh);
119     }
120   };
121   //--------------------------------------------------------------------------------
122   /*!
123    * \brief Listener of events of 3D sub-meshes computed with viscous layers.
124    * It is used to clear an inferior dim sub-meshes modified by viscous layers
125    */
126   class _ShrinkShapeListener : SMESH_subMeshEventListener
127   {
128     _ShrinkShapeListener()
129       : SMESH_subMeshEventListener(/*isDeletable=*/false,
130                                    "StdMeshers_ViscousLayers::_ShrinkShapeListener") {}
131   public:
132     static SMESH_subMeshEventListener* Get() { static _ShrinkShapeListener l; return &l; }
133     virtual void ProcessEvent(const int                       event,
134                               const int                       eventType,
135                               SMESH_subMesh*                  solidSM,
136                               SMESH_subMeshEventListenerData* data,
137                               const SMESH_Hypothesis*         hyp)
138     {
139       if ( SMESH_subMesh::COMPUTE_EVENT == eventType && solidSM->IsEmpty() && data )
140       {
141         SMESH_subMeshEventListener::ProcessEvent(event,eventType,solidSM,data,hyp);
142       }
143     }
144   };
145   //--------------------------------------------------------------------------------
146   /*!
147    * \brief Listener of events of 3D sub-meshes computed with viscous layers.
148    * It is used to store data computed by _ViscousBuilder for a sub-mesh and to
149    * delete the data as soon as it has been used
150    */
151   class _ViscousListener : SMESH_subMeshEventListener
152   {
153     _ViscousListener():
154       SMESH_subMeshEventListener(/*isDeletable=*/false,
155                                  "StdMeshers_ViscousLayers::_ViscousListener") {}
156     static SMESH_subMeshEventListener* Get() { static _ViscousListener l; return &l; }
157   public:
158     virtual void ProcessEvent(const int                       event,
159                               const int                       eventType,
160                               SMESH_subMesh*                  subMesh,
161                               SMESH_subMeshEventListenerData* data,
162                               const SMESH_Hypothesis*         hyp)
163     {
164       if ( SMESH_subMesh::COMPUTE_EVENT == eventType )
165       {
166         // delete SMESH_ProxyMesh containing temporary faces
167         subMesh->DeleteEventListener( this );
168       }
169     }
170     // Finds or creates proxy mesh of the solid
171     static _MeshOfSolid* GetSolidMesh(SMESH_Mesh*         mesh,
172                                       const TopoDS_Shape& solid,
173                                       bool                toCreate=false)
174     {
175       if ( !mesh ) return 0;
176       SMESH_subMesh* sm = mesh->GetSubMesh(solid);
177       _MeshOfSolid* data = (_MeshOfSolid*) sm->GetEventListenerData( Get() );
178       if ( !data && toCreate )
179       {
180         data = new _MeshOfSolid(mesh);
181         data->mySubMeshes.push_back( sm ); // to find SOLID by _MeshOfSolid
182         sm->SetEventListener( Get(), data, sm );
183       }
184       return data;
185     }
186     // Removes proxy mesh of the solid
187     static void RemoveSolidMesh(SMESH_Mesh* mesh, const TopoDS_Shape& solid)
188     {
189       mesh->GetSubMesh(solid)->DeleteEventListener( _ViscousListener::Get() );
190     }
191   };
192   
193   //================================================================================
194   /*!
195    * \brief sets a sub-mesh event listener to clear sub-meshes of sub-shapes of
196    * the main shape when sub-mesh of the main shape is cleared,
197    * for example to clear sub-meshes of FACEs when sub-mesh of a SOLID
198    * is cleared
199    */
200   //================================================================================
201
202   void ToClearSubWithMain( SMESH_subMesh* sub, const TopoDS_Shape& main)
203   {
204     SMESH_subMesh* mainSM = sub->GetFather()->GetSubMesh( main );
205     SMESH_subMeshEventListenerData* data =
206       mainSM->GetEventListenerData( _ShrinkShapeListener::Get());
207     if ( data )
208     {
209       if ( find( data->mySubMeshes.begin(), data->mySubMeshes.end(), sub ) ==
210            data->mySubMeshes.end())
211         data->mySubMeshes.push_back( sub );
212     }
213     else
214     {
215       data = SMESH_subMeshEventListenerData::MakeData( /*dependent=*/sub );
216       sub->SetEventListener( _ShrinkShapeListener::Get(), data, /*whereToListenTo=*/mainSM );
217     }
218   }
219   //--------------------------------------------------------------------------------
220   /*!
221    * \brief Simplex (triangle or tetrahedron) based on 1 (tria) or 2 (tet) nodes of
222    * _LayerEdge and 2 nodes of the mesh surface beening smoothed.
223    * The class is used to check validity of face or volumes around a smoothed node;
224    * it stores only 2 nodes as the other nodes are stored by _LayerEdge.
225    */
226   struct _Simplex
227   {
228     const SMDS_MeshNode *_nPrev, *_nNext; // nodes on a smoothed mesh surface
229     const SMDS_MeshNode *_nOpp; // in 2D case, a node opposite to a smoothed node in QUAD
230     _Simplex(const SMDS_MeshNode* nPrev=0,
231              const SMDS_MeshNode* nNext=0,
232              const SMDS_MeshNode* nOpp=0)
233       : _nPrev(nPrev), _nNext(nNext), _nOpp(nOpp) {}
234     bool IsForward(const SMDS_MeshNode* nSrc, const gp_XYZ* pntTgt) const
235     {
236       const double M[3][3] =
237         {{ _nNext->X() - nSrc->X(), _nNext->Y() - nSrc->Y(), _nNext->Z() - nSrc->Z() },
238          { pntTgt->X() - nSrc->X(), pntTgt->Y() - nSrc->Y(), pntTgt->Z() - nSrc->Z() },
239          { _nPrev->X() - nSrc->X(), _nPrev->Y() - nSrc->Y(), _nPrev->Z() - nSrc->Z() }};
240       double determinant = ( + M[0][0]*M[1][1]*M[2][2]
241                              + M[0][1]*M[1][2]*M[2][0]
242                              + M[0][2]*M[1][0]*M[2][1]
243                              - M[0][0]*M[1][2]*M[2][1]
244                              - M[0][1]*M[1][0]*M[2][2]
245                              - M[0][2]*M[1][1]*M[2][0]);
246       return determinant > 1e-100;
247     }
248     bool IsForward(const gp_XY&         tgtUV,
249                    const SMDS_MeshNode* smoothedNode,
250                    const TopoDS_Face&   face,
251                    SMESH_MesherHelper&  helper,
252                    const double         refSign) const
253     {
254       gp_XY prevUV = helper.GetNodeUV( face, _nPrev, smoothedNode );
255       gp_XY nextUV = helper.GetNodeUV( face, _nNext, smoothedNode );
256       gp_Vec2d v1( tgtUV, prevUV ), v2( tgtUV, nextUV );
257       double d = v1 ^ v2;
258       return d*refSign > 1e-100;
259     }
260     bool IsNeighbour(const _Simplex& other) const
261     {
262       return _nPrev == other._nNext || _nNext == other._nPrev;
263     }
264   };
265   //--------------------------------------------------------------------------------
266   /*!
267    * Structure used to take into account surface curvature while smoothing
268    */
269   struct _Curvature
270   {
271     double _r; // radius
272     double _k; // factor to correct node smoothed position
273     double _h2lenRatio; // avgNormProj / (2*avgDist)
274   public:
275     static _Curvature* New( double avgNormProj, double avgDist )
276     {
277       _Curvature* c = 0;
278       if ( fabs( avgNormProj / avgDist ) > 1./200 )
279       {
280         c = new _Curvature;
281         c->_r = avgDist * avgDist / avgNormProj;
282         c->_k = avgDist * avgDist / c->_r / c->_r;
283         c->_k *= ( c->_r < 0 ? 1/1.1 : 1.1 ); // not to be too restrictive
284         c->_h2lenRatio = avgNormProj / ( avgDist + avgDist );
285       }
286       return c;
287     }
288     double lenDelta(double len) const { return _k * ( _r + len ); }
289     double lenDeltaByDist(double dist) const { return dist * _h2lenRatio; }
290   };
291   struct _LayerEdge;
292   //--------------------------------------------------------------------------------
293   /*!
294    * Structure used to smooth a _LayerEdge (master) based on an EDGE.
295    */
296   struct _2NearEdges
297   {
298     // target nodes of 2 neighbour _LayerEdge's based on the same EDGE
299     const SMDS_MeshNode* _nodes[2];
300     // vectors from source nodes of 2 _LayerEdge's to the source node of master _LayerEdge
301     //gp_XYZ               _vec[2];
302     double               _wgt[2]; // weights of _nodes
303     _LayerEdge*          _edges[2];
304
305      // normal to plane passing through _LayerEdge._normal and tangent of EDGE
306     gp_XYZ*              _plnNorm;
307
308     _2NearEdges() { _nodes[0]=_nodes[1]=0; _plnNorm = 0; }
309     void reverse() {
310       std::swap( _nodes[0], _nodes[1] );
311       std::swap( _wgt[0], _wgt[1] );
312     }
313   };
314   //--------------------------------------------------------------------------------
315   /*!
316    * \brief Edge normal to surface, connecting a node on solid surface (_nodes[0])
317    * and a node of the most internal layer (_nodes.back())
318    */
319   struct _LayerEdge
320   {
321     vector< const SMDS_MeshNode*> _nodes;
322
323     gp_XYZ              _normal; // to solid surface
324     vector<gp_XYZ>      _pos; // points computed during inflation
325     double              _len; // length achived with the last step
326     double              _cosin; // of angle (_normal ^ surface)
327     double              _lenFactor; // to compute _len taking _cosin into account
328
329     // face or edge w/o layer along or near which _LayerEdge is inflated
330     TopoDS_Shape        _sWOL;
331     // simplices connected to the source node (_nodes[0]);
332     // used for smoothing and quality check of _LayerEdge's based on the FACE
333     vector<_Simplex>    _simplices;
334     // data for smoothing of _LayerEdge's based on the EDGE
335     _2NearEdges*        _2neibors;
336
337     _Curvature*         _curvature;
338     // TODO:: detele _Curvature, _plnNorm
339
340     void SetNewLength( double len, SMESH_MesherHelper& helper );
341     bool SetNewLength2d( Handle(Geom_Surface)& surface,
342                          const TopoDS_Face&    F,
343                          SMESH_MesherHelper&   helper );
344     void SetDataByNeighbors( const SMDS_MeshNode* n1,
345                              const SMDS_MeshNode* n2,
346                              SMESH_MesherHelper&  helper);
347     void InvalidateStep( int curStep );
348     bool Smooth(int& badNb);
349     bool SmoothOnEdge(Handle(Geom_Surface)& surface,
350                       const TopoDS_Face&    F,
351                       SMESH_MesherHelper&   helper);
352     bool FindIntersection( SMESH_ElementSearcher&   searcher,
353                            double &                 distance,
354                            const double&            epsilon,
355                            const SMDS_MeshElement** face = 0);
356     bool SegTriaInter( const gp_Ax1&        lastSegment,
357                        const SMDS_MeshNode* n0,
358                        const SMDS_MeshNode* n1,
359                        const SMDS_MeshNode* n2,
360                        double&              dist,
361                        const double&        epsilon) const;
362     gp_Ax1 LastSegment(double& segLen) const;
363     bool IsOnEdge() const { return _2neibors; }
364     void Copy( _LayerEdge& other, SMESH_MesherHelper& helper );
365     void SetCosin( double cosin );
366   };
367   struct _LayerEdgeCmp
368   {
369     bool operator () (const _LayerEdge* e1, const _LayerEdge* e2) const
370     {
371       const bool cmpNodes = ( e1 && e2 && e1->_nodes.size() && e2->_nodes.size() );
372       return cmpNodes ? ( e1->_nodes[0]->GetID() < e2->_nodes[0]->GetID()) : ( e1 < e2 );
373     }
374   };
375   //--------------------------------------------------------------------------------
376
377   typedef map< const SMDS_MeshNode*, _LayerEdge*, TIDCompare > TNode2Edge;
378   
379   //--------------------------------------------------------------------------------
380   /*!
381    * \brief Data of a SOLID
382    */
383   struct _SolidData
384   {
385     TopoDS_Shape                    _solid;
386     const StdMeshers_ViscousLayers* _hyp;
387     TopoDS_Shape                    _hypShape;
388     _MeshOfSolid*                   _proxyMesh;
389     set<TGeomID>                    _reversedFaceIds;
390     set<TGeomID>                    _ignoreFaceIds;
391
392     double                          _stepSize, _stepSizeCoeff;
393     const SMDS_MeshNode*            _stepSizeNodes[2];
394
395     TNode2Edge                      _n2eMap;
396     // edges of _n2eMap. We keep same data in two containers because
397     // iteration over the map is 5 time longer than over the vector
398     vector< _LayerEdge* >           _edges;
399
400     // key:   an id of shape (EDGE or VERTEX) shared by a FACE with
401     //        layers and a FACE w/o layers
402     // value: the shape (FACE or EDGE) to shrink mesh on.
403     //       _LayerEdge's basing on nodes on key shape are inflated along the value shape
404     map< TGeomID, TopoDS_Shape >     _shrinkShape2Shape;
405
406     // FACE's WOL, srink on which is forbiden due to algo on the adjacent SOLID
407     set< TGeomID >                   _noShrinkFaces;
408
409     // <EDGE to smooth on> to <it's curve>
410     map< TGeomID,Handle(Geom_Curve)> _edge2curve;
411
412     // end indices in _edges of _LayerEdge on one shape to smooth
413     vector< int >                    _endEdgeToSmooth;
414
415     double                           _epsilon; // precision for SegTriaInter()
416
417     int                              _index; // for debug
418
419     _SolidData(const TopoDS_Shape&             s=TopoDS_Shape(),
420                const StdMeshers_ViscousLayers* h=0,
421                const TopoDS_Shape&             hs=TopoDS_Shape(),
422                _MeshOfSolid*                   m=0)
423       :_solid(s), _hyp(h), _hypShape(hs), _proxyMesh(m) {}
424     ~_SolidData();
425
426     Handle(Geom_Curve) CurveForSmooth( const TopoDS_Edge&    E,
427                                        const int             iFrom,
428                                        const int             iTo,
429                                        Handle(Geom_Surface)& surface,
430                                        const TopoDS_Face&    F,
431                                        SMESH_MesherHelper&   helper);
432   };
433   //--------------------------------------------------------------------------------
434   /*!
435    * \brief Data of node on a shrinked FACE
436    */
437   struct _SmoothNode
438   {
439     const SMDS_MeshNode*         _node;
440     //vector<const SMDS_MeshNode*> _nodesAround;
441     vector<_Simplex>             _simplices; // for quality check
442
443     enum SmoothType { LAPLACIAN, CENTROIDAL, ANGULAR, TFI };
444
445     bool Smooth(int&                  badNb,
446                 Handle(Geom_Surface)& surface,
447                 SMESH_MesherHelper&   helper,
448                 const double          refSign,
449                 SmoothType            how,
450                 bool                  set3D);
451
452     gp_XY computeAngularPos(vector<gp_XY>& uv,
453                             const gp_XY&   uvToFix,
454                             const double   refSign );
455   };
456   //--------------------------------------------------------------------------------
457   /*!
458    * \brief Builder of viscous layers
459    */
460   class _ViscousBuilder
461   {
462   public:
463     _ViscousBuilder();
464     // does it's job
465     SMESH_ComputeErrorPtr Compute(SMESH_Mesh&         mesh,
466                                   const TopoDS_Shape& shape);
467
468     // restore event listeners used to clear an inferior dim sub-mesh modified by viscous layers
469     void RestoreListeners();
470
471     // computes SMESH_ProxyMesh::SubMesh::_n2n;
472     bool MakeN2NMap( _MeshOfSolid* pm );
473
474   private:
475
476     bool findSolidsWithLayers();
477     bool findFacesWithLayers();
478     bool makeLayer(_SolidData& data);
479     bool setEdgeData(_LayerEdge& edge, const set<TGeomID>& subIds,
480                      SMESH_MesherHelper& helper, _SolidData& data);
481     bool findNeiborsOnEdge(const _LayerEdge*     edge,
482                            const SMDS_MeshNode*& n1,
483                            const SMDS_MeshNode*& n2,
484                            _SolidData&           data);
485     void getSimplices( const SMDS_MeshNode* node, vector<_Simplex>& simplices,
486                        const set<TGeomID>& ingnoreShapes,
487                        const _SolidData*   dataToCheckOri = 0,
488                        const bool          toSort = false);
489     bool sortEdges( _SolidData&                    data,
490                     vector< vector<_LayerEdge*> >& edgesByGeom);
491     void limitStepSize( _SolidData&             data,
492                         const SMDS_MeshElement* face,
493                         const double            cosin);
494     void limitStepSize( _SolidData& data, const double minSize);
495     bool inflate(_SolidData& data);
496     bool smoothAndCheck(_SolidData& data, const int nbSteps, double & distToIntersection);
497     bool smoothAnalyticEdge( _SolidData&           data,
498                              const int             iFrom,
499                              const int             iTo,
500                              Handle(Geom_Surface)& surface,
501                              const TopoDS_Face&    F,
502                              SMESH_MesherHelper&   helper);
503     bool updateNormals( _SolidData& data, SMESH_MesherHelper& helper );
504     bool refine(_SolidData& data);
505     bool shrink();
506     bool prepareEdgeToShrink( _LayerEdge& edge, const TopoDS_Face& F,
507                               SMESH_MesherHelper& helper,
508                               const SMESHDS_SubMesh* faceSubMesh );
509     void fixBadFaces(const TopoDS_Face&          F,
510                      SMESH_MesherHelper&         helper,
511                      const bool                  is2D,
512                      const int                   step,
513                      set<const SMDS_MeshNode*> * involvedNodes=NULL);
514     bool addBoundaryElements();
515
516     bool error( const string& text, int solidID=-1 );
517     SMESHDS_Mesh* getMeshDS() { return _mesh->GetMeshDS(); }
518
519     // debug
520     void makeGroupOfLE();
521
522     SMESH_Mesh*           _mesh;
523     SMESH_ComputeErrorPtr _error;
524
525     vector< _SolidData >  _sdVec;
526     int                   _tmpFaceID;
527   };
528   //--------------------------------------------------------------------------------
529   /*!
530    * \brief Shrinker of nodes on the EDGE
531    */
532   class _Shrinker1D
533   {
534     vector<double>                _initU;
535     vector<double>                _normPar;
536     vector<const SMDS_MeshNode*>  _nodes;
537     const _LayerEdge*             _edges[2];
538     bool                          _done;
539   public:
540     void AddEdge( const _LayerEdge* e, SMESH_MesherHelper& helper );
541     void Compute(bool set3D, SMESH_MesherHelper& helper);
542     void RestoreParams();
543     void SwapSrcTgtNodes(SMESHDS_Mesh* mesh);
544   };
545   //--------------------------------------------------------------------------------
546   /*!
547    * \brief Class of temporary mesh face.
548    * We can't use SMDS_FaceOfNodes since it's impossible to set it's ID which is
549    * needed because SMESH_ElementSearcher internaly uses set of elements sorted by ID
550    */
551   struct TmpMeshFace : public SMDS_MeshElement
552   {
553     vector<const SMDS_MeshNode* > _nn;
554     TmpMeshFace( const vector<const SMDS_MeshNode*>& nodes, int id):
555       SMDS_MeshElement(id), _nn(nodes) {}
556     virtual const SMDS_MeshNode* GetNode(const int ind) const { return _nn[ind]; }
557     virtual SMDSAbs_ElementType  GetType() const              { return SMDSAbs_Face; }
558     virtual vtkIdType GetVtkType() const                      { return -1; }
559     virtual SMDSAbs_EntityType   GetEntityType() const        { return SMDSEntity_Last; }
560     virtual SMDSAbs_GeometryType GetGeomType() const          { return SMDSGeom_TRIANGLE; }
561     virtual SMDS_ElemIteratorPtr elementsIterator(SMDSAbs_ElementType) const
562     { return SMDS_ElemIteratorPtr( new SMDS_NodeVectorElemIterator( _nn.begin(), _nn.end()));}
563   };
564   //--------------------------------------------------------------------------------
565   /*!
566    * \brief Class of temporary mesh face storing _LayerEdge it's based on
567    */
568   struct TmpMeshFaceOnEdge : public TmpMeshFace
569   {
570     _LayerEdge *_le1, *_le2;
571     TmpMeshFaceOnEdge( _LayerEdge* le1, _LayerEdge* le2, int ID ):
572       TmpMeshFace( vector<const SMDS_MeshNode*>(4), ID ), _le1(le1), _le2(le2)
573     {
574       _nn[0]=_le1->_nodes[0];
575       _nn[1]=_le1->_nodes.back();
576       _nn[2]=_le2->_nodes.back();
577       _nn[3]=_le2->_nodes[0];
578     }
579   };
580   //--------------------------------------------------------------------------------
581   /*!
582    * \brief Retriever of node coordinates either directly of from a surface by node UV.
583    * \warning Location of a surface is ignored
584    */
585   struct NodeCoordHelper
586   {
587     SMESH_MesherHelper&        _helper;
588     const TopoDS_Face&         _face;
589     Handle(Geom_Surface)       _surface;
590     gp_XYZ (NodeCoordHelper::* _fun)(const SMDS_MeshNode* n) const;
591
592     NodeCoordHelper(const TopoDS_Face& F, SMESH_MesherHelper& helper, bool is2D)
593       : _helper( helper ), _face( F )
594     {
595       if ( is2D )
596       {
597         TopLoc_Location loc;
598         _surface = BRep_Tool::Surface( _face, loc );
599       }
600       if ( _surface.IsNull() )
601         _fun = & NodeCoordHelper::direct;
602       else
603         _fun = & NodeCoordHelper::byUV;
604     }
605     gp_XYZ operator()(const SMDS_MeshNode* n) const { return (this->*_fun)( n ); }
606
607   private:
608     gp_XYZ direct(const SMDS_MeshNode* n) const
609     {
610       return SMESH_TNodeXYZ( n );
611     }
612     gp_XYZ byUV  (const SMDS_MeshNode* n) const
613     {
614       gp_XY uv = _helper.GetNodeUV( _face, n );
615       return _surface->Value( uv.X(), uv.Y() ).XYZ();
616     }
617   };
618 } // namespace VISCOUS_3D
619
620 //================================================================================
621 // StdMeshers_ViscousLayers hypothesis
622 //
623 StdMeshers_ViscousLayers::StdMeshers_ViscousLayers(int hypId, int studyId, SMESH_Gen* gen)
624   :SMESH_Hypothesis(hypId, studyId, gen),
625    _isToIgnoreShapes(1), _nbLayers(1), _thickness(1), _stretchFactor(1)
626 {
627   _name = StdMeshers_ViscousLayers::GetHypType();
628   _param_algo_dim = -3; // auxiliary hyp used by 3D algos
629 } // --------------------------------------------------------------------------------
630 void StdMeshers_ViscousLayers::SetBndShapes(const std::vector<int>& faceIds, bool toIgnore)
631 {
632   if ( faceIds != _shapeIds )
633     _shapeIds = faceIds, NotifySubMeshesHypothesisModification();
634   if ( _isToIgnoreShapes != toIgnore )
635     _isToIgnoreShapes = toIgnore, NotifySubMeshesHypothesisModification();
636 } // --------------------------------------------------------------------------------
637 void StdMeshers_ViscousLayers::SetTotalThickness(double thickness)
638 {
639   if ( thickness != _thickness )
640     _thickness = thickness, NotifySubMeshesHypothesisModification();
641 } // --------------------------------------------------------------------------------
642 void StdMeshers_ViscousLayers::SetNumberLayers(int nb)
643 {
644   if ( _nbLayers != nb )
645     _nbLayers = nb, NotifySubMeshesHypothesisModification();
646 } // --------------------------------------------------------------------------------
647 void StdMeshers_ViscousLayers::SetStretchFactor(double factor)
648 {
649   if ( _stretchFactor != factor )
650     _stretchFactor = factor, NotifySubMeshesHypothesisModification();
651 } // --------------------------------------------------------------------------------
652 SMESH_ProxyMesh::Ptr
653 StdMeshers_ViscousLayers::Compute(SMESH_Mesh&         theMesh,
654                                   const TopoDS_Shape& theShape,
655                                   const bool          toMakeN2NMap) const
656 {
657   using namespace VISCOUS_3D;
658   _ViscousBuilder bulder;
659   SMESH_ComputeErrorPtr err = bulder.Compute( theMesh, theShape );
660   if ( err && !err->IsOK() )
661     return SMESH_ProxyMesh::Ptr();
662
663   vector<SMESH_ProxyMesh::Ptr> components;
664   TopExp_Explorer exp( theShape, TopAbs_SOLID );
665   for ( ; exp.More(); exp.Next() )
666   {
667     if ( _MeshOfSolid* pm =
668          _ViscousListener::GetSolidMesh( &theMesh, exp.Current(), /*toCreate=*/false))
669     {
670       if ( toMakeN2NMap && !pm->_n2nMapComputed )
671         if ( !bulder.MakeN2NMap( pm ))
672           return SMESH_ProxyMesh::Ptr();
673       components.push_back( SMESH_ProxyMesh::Ptr( pm ));
674       pm->myIsDeletable = false; // it will de deleted by boost::shared_ptr
675     }
676     _ViscousListener::RemoveSolidMesh ( &theMesh, exp.Current() );
677   }
678   switch ( components.size() )
679   {
680   case 0: break;
681
682   case 1: return components[0];
683
684   default: return SMESH_ProxyMesh::Ptr( new SMESH_ProxyMesh( components ));
685   }
686   return SMESH_ProxyMesh::Ptr();
687 } // --------------------------------------------------------------------------------
688 std::ostream & StdMeshers_ViscousLayers::SaveTo(std::ostream & save)
689 {
690   save << " " << _nbLayers
691        << " " << _thickness
692        << " " << _stretchFactor
693        << " " << _shapeIds.size();
694   for ( size_t i = 0; i < _shapeIds.size(); ++i )
695     save << " " << _shapeIds[i];
696   save << " " << !_isToIgnoreShapes; // negate to keep the behavior in old studies.
697   return save;
698 } // --------------------------------------------------------------------------------
699 std::istream & StdMeshers_ViscousLayers::LoadFrom(std::istream & load)
700 {
701   int nbFaces, faceID, shapeToTreat;
702   load >> _nbLayers >> _thickness >> _stretchFactor >> nbFaces;
703   while ( _shapeIds.size() < nbFaces && load >> faceID )
704     _shapeIds.push_back( faceID );
705   if ( load >> shapeToTreat )
706     _isToIgnoreShapes = !shapeToTreat;
707   else
708     _isToIgnoreShapes = true; // old behavior
709   return load;
710 } // --------------------------------------------------------------------------------
711 bool StdMeshers_ViscousLayers::SetParametersByMesh(const SMESH_Mesh*   theMesh,
712                                                    const TopoDS_Shape& theShape)
713 {
714   // TODO
715   return false;
716 }
717 // END StdMeshers_ViscousLayers hypothesis
718 //================================================================================
719
720 namespace
721 {
722   gp_XYZ getEdgeDir( const TopoDS_Edge& E, const TopoDS_Vertex& fromV )
723   {
724     gp_Vec dir;
725     double f,l;
726     Handle(Geom_Curve) c = BRep_Tool::Curve( E, f, l );
727     gp_Pnt p = BRep_Tool::Pnt( fromV );
728     double distF = p.SquareDistance( c->Value( f ));
729     double distL = p.SquareDistance( c->Value( l ));
730     c->D1(( distF < distL ? f : l), p, dir );
731     if ( distL < distF ) dir.Reverse();
732     return dir.XYZ();
733   }
734   //--------------------------------------------------------------------------------
735   gp_XYZ getEdgeDir( const TopoDS_Edge& E, const SMDS_MeshNode* atNode,
736                      SMESH_MesherHelper& helper)
737   {
738     gp_Vec dir;
739     double f,l; gp_Pnt p;
740     Handle(Geom_Curve) c = BRep_Tool::Curve( E, f, l );
741     double u = helper.GetNodeU( E, atNode );
742     c->D1( u, p, dir );
743     return dir.XYZ();
744   }
745   //--------------------------------------------------------------------------------
746   gp_XYZ getFaceDir( const TopoDS_Face& F, const TopoDS_Edge& fromE,
747                      const SMDS_MeshNode* node, SMESH_MesherHelper& helper, bool& ok)
748   {
749     gp_XY uv = helper.GetNodeUV( F, node, 0, &ok );
750     Handle(Geom_Surface) surface = BRep_Tool::Surface( F );
751     gp_Pnt p; gp_Vec du, dv, norm;
752     surface->D1( uv.X(),uv.Y(), p, du,dv );
753     norm = du ^ dv;
754
755     double f,l;
756     Handle(Geom_Curve) c = BRep_Tool::Curve( fromE, f, l );
757     double u = helper.GetNodeU( fromE, node, 0, &ok );
758     c->D1( u, p, du );
759     TopAbs_Orientation o = helper.GetSubShapeOri( F.Oriented(TopAbs_FORWARD), fromE);
760     if ( o == TopAbs_REVERSED )
761       du.Reverse();
762
763     gp_Vec dir = norm ^ du;
764
765     if ( node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX &&
766          helper.IsClosedEdge( fromE ))
767     {
768       if ( fabs(u-f) < fabs(u-l )) c->D1( l, p, dv );
769       else                         c->D1( f, p, dv );
770       if ( o == TopAbs_REVERSED )
771         dv.Reverse();
772       gp_Vec dir2 = norm ^ dv;
773       dir = dir.Normalized() + dir2.Normalized();
774     }
775     return dir.XYZ();
776   }
777   //--------------------------------------------------------------------------------
778   gp_XYZ getFaceDir( const TopoDS_Face& F, const TopoDS_Vertex& fromV,
779                      const SMDS_MeshNode* node, SMESH_MesherHelper& helper,
780                      bool& ok, double* cosin=0)
781   {
782     double f,l; TopLoc_Location loc;
783     vector< TopoDS_Edge > edges; // sharing a vertex
784     PShapeIteratorPtr eIt = helper.GetAncestors( fromV, *helper.GetMesh(), TopAbs_EDGE);
785     while ( eIt->more())
786     {
787       const TopoDS_Edge* e = static_cast<const TopoDS_Edge*>( eIt->next() );
788       if ( helper.IsSubShape( *e, F ) && !BRep_Tool::Curve( *e, loc,f,l).IsNull() )
789         edges.push_back( *e );
790     }
791     gp_XYZ dir(0,0,0);
792     if ( !( ok = ( edges.size() > 0 ))) return dir;
793     // get average dir of edges going fromV
794     gp_XYZ edgeDir;
795     //if ( edges.size() > 1 )
796       for ( size_t i = 0; i < edges.size(); ++i )
797       {
798         edgeDir = getEdgeDir( edges[i], fromV );
799         double size2 = edgeDir.SquareModulus();
800         if ( size2 > numeric_limits<double>::min() )
801           edgeDir /= sqrt( size2 );
802         else
803           ok = false;
804         dir += edgeDir;
805       }
806     gp_XYZ fromEdgeDir = getFaceDir( F, edges[0], node, helper, ok );
807     if ( edges.size() == 1 )
808       dir = fromEdgeDir;
809     else if ( dir.SquareModulus() < 0.1 ) // ~< 20 degrees
810       dir = fromEdgeDir + getFaceDir( F, edges[1], node, helper, ok );
811     else if ( dir * fromEdgeDir < 0 )
812       dir *= -1;
813     if ( ok )
814     {
815       //dir /= edges.size();
816       if ( cosin ) {
817         double angle = gp_Vec( edgeDir ).Angle( dir );
818         *cosin = cos( angle );
819       }
820     }
821     return dir;
822   }
823   //================================================================================
824   /*!
825    * \brief Returns true if a FACE is bound by a concave EDGE
826    */
827   //================================================================================
828
829   bool isConcave( const TopoDS_Face& F, SMESH_MesherHelper& helper )
830   {
831     // if ( helper.Count( F, TopAbs_WIRE, /*useMap=*/false) > 1 )
832     //   return true;
833     gp_Vec2d drv1, drv2;
834     gp_Pnt2d p;
835     TopExp_Explorer eExp( F.Oriented( TopAbs_FORWARD ), TopAbs_EDGE );
836     for ( ; eExp.More(); eExp.Next() )
837     {
838       const TopoDS_Edge& E = TopoDS::Edge( eExp.Current() );
839       if ( SMESH_Algo::isDegenerated( E )) continue;
840       // check if 2D curve is concave
841       BRepAdaptor_Curve2d curve( E, F );
842       const int nbIntervals = curve.NbIntervals( GeomAbs_C2 );
843       TColStd_Array1OfReal intervals(1, nbIntervals + 1 );
844       curve.Intervals( intervals, GeomAbs_C2 );
845       bool isConvex = true;
846       for ( int i = 1; i <= nbIntervals && isConvex; ++i )
847       {
848         double u1 = intervals( i );
849         double u2 = intervals( i+1 );
850         curve.D2( 0.5*( u1+u2 ), p, drv1, drv2 );
851         double cross = drv2 ^ drv1;
852         if ( E.Orientation() == TopAbs_REVERSED )
853           cross = -cross;
854         isConvex = ( cross > 0.1 ); //-1e-9 );
855       }
856       // check if concavity is strong enough to care about it
857       //const double maxAngle = 5 * Standard_PI180;
858       if ( !isConvex )
859       {
860         //cout << "Concave FACE " << helper.GetMeshDS()->ShapeToIndex( F ) << endl;
861         return true;
862         // map< double, const SMDS_MeshNode* > u2nodes;
863         // if ( !SMESH_Algo::GetSortedNodesOnEdge( helper.GetMeshDS(), E,
864         //                                         /*ignoreMedium=*/true, u2nodes))
865         //   continue;
866         // map< double, const SMDS_MeshNode* >::iterator u2n = u2nodes.begin();
867         // gp_Pnt2d uvPrev = helper.GetNodeUV( F, u2n->second );
868         // double    uPrev = u2n->first;
869         // for ( ++u2n; u2n != u2nodes.end(); ++u2n )
870         // {
871         //   gp_Pnt2d uv = helper.GetNodeUV( F, u2n->second );
872         //   gp_Vec2d segmentDir( uvPrev, uv );
873         //   curve.D1( uPrev, p, drv1 );
874         //   try {
875         //     if ( fabs( segmentDir.Angle( drv1 )) > maxAngle )
876         //       return true;
877         //   }
878         //   catch ( ... ) {}
879         //   uvPrev = uv;
880         //   uPrev = u2n->first;
881         // }
882       }
883     }
884     // check angles at VERTEXes
885     TError error;
886     TSideVector wires = StdMeshers_FaceSide::GetFaceWires( F, *helper.GetMesh(), 0, error );
887     for ( size_t iW = 0; iW < wires.size(); ++iW )
888     {
889       const int nbEdges = wires[iW]->NbEdges();
890       if ( nbEdges < 2 && SMESH_Algo::isDegenerated( wires[iW]->Edge(0)))
891         continue;
892       for ( int iE1 = 0; iE1 < nbEdges; ++iE1 )
893       {
894         if ( SMESH_Algo::isDegenerated( wires[iW]->Edge( iE1 ))) continue;
895         int iE2 = ( iE1 + 1 ) % nbEdges;
896         while ( SMESH_Algo::isDegenerated( wires[iW]->Edge( iE2 )))
897           iE2 = ( iE2 + 1 ) % nbEdges;
898         double angle = helper.GetAngle( wires[iW]->Edge( iE1 ),
899                                         wires[iW]->Edge( iE2 ), F );
900         if ( angle < -5. * M_PI / 180. )
901           return true;
902       }
903     }
904     return false;
905   }
906   //--------------------------------------------------------------------------------
907   // DEBUG. Dump intermediate node positions into a python script
908 #ifdef __myDEBUG
909   ofstream* py;
910   struct PyDump {
911     PyDump() {
912       const char* fname = "/tmp/viscous.py";
913       cout << "execfile('"<<fname<<"')"<<endl;
914       py = new ofstream(fname);
915       *py << "import SMESH" << endl
916           << "from salome.smesh import smeshBuilder" << endl
917           << "smesh  = smeshBuilder.New(salome.myStudy)" << endl
918           << "meshSO = smesh.GetCurrentStudy().FindObjectID('0:1:2:3')" << endl
919           << "mesh   = smesh.Mesh( meshSO.GetObject() )"<<endl;
920     }
921     void Finish() {
922       if (py)
923         *py << "mesh.MakeGroup('Viscous Prisms',SMESH.VOLUME,SMESH.FT_ElemGeomType,'=',SMESH.Geom_PENTA)"<<endl;
924       delete py; py=0;
925     }
926     ~PyDump() { Finish(); }
927   };
928 #define dumpFunction(f) { _dumpFunction(f, __LINE__);}
929 #define dumpMove(n)     { _dumpMove(n, __LINE__);}
930 #define dumpCmd(txt)    { _dumpCmd(txt, __LINE__);}
931   void _dumpFunction(const string& fun, int ln)
932   { if (py) *py<< "def "<<fun<<"(): # "<< ln <<endl; cout<<fun<<"()"<<endl;}
933   void _dumpMove(const SMDS_MeshNode* n, int ln)
934   { if (py) *py<< "  mesh.MoveNode( "<<n->GetID()<< ", "<< n->X()
935                << ", "<<n->Y()<<", "<< n->Z()<< ")\t\t # "<< ln <<endl; }
936   void _dumpCmd(const string& txt, int ln)
937   { if (py) *py<< "  "<<txt<<" # "<< ln <<endl; }
938   void dumpFunctionEnd()
939   { if (py) *py<< "  return"<< endl; }
940   void dumpChangeNodes( const SMDS_MeshElement* f )
941   { if (py) { *py<< "  mesh.ChangeElemNodes( " << f->GetID()<<", [";
942       for ( int i=1; i < f->NbNodes(); ++i ) *py << f->GetNode(i-1)->GetID()<<", ";
943       *py << f->GetNode( f->NbNodes()-1 )->GetID() << " ])"<< endl; }}
944 #else
945   struct PyDump { void Finish() {} };
946 #define dumpFunction(f) f
947 #define dumpMove(n)
948 #define dumpCmd(txt)
949 #define dumpFunctionEnd()
950 #define dumpChangeNodes(f)
951 #endif
952 }
953
954 using namespace VISCOUS_3D;
955
956 //================================================================================
957 /*!
958  * \brief Constructor of _ViscousBuilder
959  */
960 //================================================================================
961
962 _ViscousBuilder::_ViscousBuilder()
963 {
964   _error = SMESH_ComputeError::New(COMPERR_OK);
965   _tmpFaceID = 0;
966 }
967
968 //================================================================================
969 /*!
970  * \brief Stores error description and returns false
971  */
972 //================================================================================
973
974 bool _ViscousBuilder::error(const string& text, int solidId )
975 {
976   _error->myName    = COMPERR_ALGO_FAILED;
977   _error->myComment = string("Viscous layers builder: ") + text;
978   if ( _mesh )
979   {
980     SMESH_subMesh* sm = _mesh->GetSubMeshContaining( solidId );
981     if ( !sm && !_sdVec.empty() )
982       sm = _mesh->GetSubMeshContaining( _sdVec[0]._index );
983     if ( sm && sm->GetSubShape().ShapeType() == TopAbs_SOLID )
984     {
985       SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
986       if ( smError && smError->myAlgo )
987         _error->myAlgo = smError->myAlgo;
988       smError = _error;
989     }
990   }
991   makeGroupOfLE(); // debug
992
993   return false;
994 }
995
996 //================================================================================
997 /*!
998  * \brief At study restoration, restore event listeners used to clear an inferior
999  *  dim sub-mesh modified by viscous layers
1000  */
1001 //================================================================================
1002
1003 void _ViscousBuilder::RestoreListeners()
1004 {
1005   // TODO
1006 }
1007
1008 //================================================================================
1009 /*!
1010  * \brief computes SMESH_ProxyMesh::SubMesh::_n2n
1011  */
1012 //================================================================================
1013
1014 bool _ViscousBuilder::MakeN2NMap( _MeshOfSolid* pm )
1015 {
1016   SMESH_subMesh* solidSM = pm->mySubMeshes.front();
1017   TopExp_Explorer fExp( solidSM->GetSubShape(), TopAbs_FACE );
1018   for ( ; fExp.More(); fExp.Next() )
1019   {
1020     SMESHDS_SubMesh* srcSmDS = pm->GetMeshDS()->MeshElements( fExp.Current() );
1021     const SMESH_ProxyMesh::SubMesh* prxSmDS = pm->GetProxySubMesh( fExp.Current() );
1022
1023     if ( !srcSmDS || !prxSmDS || !srcSmDS->NbElements() || !prxSmDS->NbElements() )
1024       continue;
1025     if ( srcSmDS->GetElements()->next() == prxSmDS->GetElements()->next())
1026       continue;
1027
1028     if ( srcSmDS->NbElements() != prxSmDS->NbElements() )
1029       return error( "Different nb elements in a source and a proxy sub-mesh", solidSM->GetId());
1030
1031     SMDS_ElemIteratorPtr srcIt = srcSmDS->GetElements();
1032     SMDS_ElemIteratorPtr prxIt = prxSmDS->GetElements();
1033     while( prxIt->more() )
1034     {
1035       const SMDS_MeshElement* fSrc = srcIt->next();
1036       const SMDS_MeshElement* fPrx = prxIt->next();
1037       if ( fSrc->NbNodes() != fPrx->NbNodes())
1038         return error( "Different elements in a source and a proxy sub-mesh", solidSM->GetId());
1039       for ( int i = 0 ; i < fPrx->NbNodes(); ++i )
1040         pm->setNode2Node( fSrc->GetNode(i), fPrx->GetNode(i), prxSmDS );
1041     }
1042   }
1043   pm->_n2nMapComputed = true;
1044   return true;
1045 }
1046
1047 //================================================================================
1048 /*!
1049  * \brief Does its job
1050  */
1051 //================================================================================
1052
1053 SMESH_ComputeErrorPtr _ViscousBuilder::Compute(SMESH_Mesh&         theMesh,
1054                                                const TopoDS_Shape& theShape)
1055 {
1056   // TODO: set priority of solids during Gen::Compute()
1057
1058   _mesh = & theMesh;
1059
1060   // check if proxy mesh already computed
1061   TopExp_Explorer exp( theShape, TopAbs_SOLID );
1062   if ( !exp.More() )
1063     return error("No SOLID's in theShape"), _error;
1064
1065   if ( _ViscousListener::GetSolidMesh( _mesh, exp.Current(), /*toCreate=*/false))
1066     return SMESH_ComputeErrorPtr(); // everything already computed
1067
1068   PyDump debugDump;
1069
1070   // TODO: ignore already computed SOLIDs 
1071   if ( !findSolidsWithLayers())
1072     return _error;
1073
1074   if ( !findFacesWithLayers() )
1075     return _error;
1076
1077   for ( size_t i = 0; i < _sdVec.size(); ++i )
1078   {
1079     if ( ! makeLayer(_sdVec[i]) )
1080       return _error;
1081
1082     if ( _sdVec[i]._edges.size() == 0 )
1083       continue;
1084     
1085     if ( ! inflate(_sdVec[i]) )
1086       return _error;
1087
1088     if ( ! refine(_sdVec[i]) )
1089       return _error;
1090   }
1091   if ( !shrink() )
1092     return _error;
1093
1094   addBoundaryElements();
1095
1096   makeGroupOfLE(); // debug
1097   debugDump.Finish();
1098
1099   return _error;
1100 }
1101
1102 //================================================================================
1103 /*!
1104  * \brief Finds SOLIDs to compute using viscous layers. Fills _sdVec
1105  */
1106 //================================================================================
1107
1108 bool _ViscousBuilder::findSolidsWithLayers()
1109 {
1110   // get all solids
1111   TopTools_IndexedMapOfShape allSolids;
1112   TopExp::MapShapes( _mesh->GetShapeToMesh(), TopAbs_SOLID, allSolids );
1113   _sdVec.reserve( allSolids.Extent());
1114
1115   SMESH_Gen* gen = _mesh->GetGen();
1116   SMESH_HypoFilter filter;
1117   for ( int i = 1; i <= allSolids.Extent(); ++i )
1118   {
1119     // find StdMeshers_ViscousLayers hyp assigned to the i-th solid
1120     SMESH_Algo* algo = gen->GetAlgo( *_mesh, allSolids(i) );
1121     if ( !algo ) continue;
1122     // TODO: check if algo is hidden
1123     const list <const SMESHDS_Hypothesis *> & allHyps =
1124       algo->GetUsedHypothesis(*_mesh, allSolids(i), /*ignoreAuxiliary=*/false);
1125     list< const SMESHDS_Hypothesis *>::const_iterator hyp = allHyps.begin();
1126     const StdMeshers_ViscousLayers* viscHyp = 0;
1127     for ( ; hyp != allHyps.end() && !viscHyp; ++hyp )
1128       viscHyp = dynamic_cast<const StdMeshers_ViscousLayers*>( *hyp );
1129     if ( viscHyp )
1130     {
1131       TopoDS_Shape hypShape;
1132       filter.Init( filter.Is( viscHyp ));
1133       _mesh->GetHypothesis( allSolids(i), filter, true, &hypShape );
1134
1135       _MeshOfSolid* proxyMesh = _ViscousListener::GetSolidMesh( _mesh,
1136                                                                 allSolids(i),
1137                                                                 /*toCreate=*/true);
1138       _sdVec.push_back( _SolidData( allSolids(i), viscHyp, hypShape, proxyMesh ));
1139       _sdVec.back()._index = getMeshDS()->ShapeToIndex( allSolids(i));
1140     }
1141   }
1142   if ( _sdVec.empty() )
1143     return error
1144       ( SMESH_Comment(StdMeshers_ViscousLayers::GetHypType()) << " hypothesis not found",0);
1145
1146   return true;
1147 }
1148
1149 //================================================================================
1150 /*!
1151  * \brief 
1152  */
1153 //================================================================================
1154
1155 bool _ViscousBuilder::findFacesWithLayers()
1156 {
1157   SMESH_MesherHelper helper( *_mesh );
1158   TopExp_Explorer exp;
1159   TopTools_IndexedMapOfShape solids;
1160
1161   // collect all faces to ignore defined by hyp
1162   for ( size_t i = 0; i < _sdVec.size(); ++i )
1163   {
1164     solids.Add( _sdVec[i]._solid );
1165
1166     vector<TGeomID> ids = _sdVec[i]._hyp->GetBndShapes();
1167     if ( _sdVec[i]._hyp->IsToIgnoreShapes() ) // FACEs to ignore are given
1168     {
1169       for ( size_t ii = 0; ii < ids.size(); ++ii )
1170       {
1171         const TopoDS_Shape& s = getMeshDS()->IndexToShape( ids[ii] );
1172         if ( !s.IsNull() && s.ShapeType() == TopAbs_FACE )
1173           _sdVec[i]._ignoreFaceIds.insert( ids[ii] );
1174       }
1175     }
1176     else // FACEs with layers are given
1177     {
1178       exp.Init( _sdVec[i]._solid, TopAbs_FACE );
1179       for ( ; exp.More(); exp.Next() )
1180       {
1181         TGeomID faceInd = getMeshDS()->ShapeToIndex( exp.Current() );
1182         if ( find( ids.begin(), ids.end(), faceInd ) == ids.end() )
1183           _sdVec[i]._ignoreFaceIds.insert( faceInd );
1184       }
1185     }
1186
1187     // ignore internal FACEs if inlets and outlets are specified
1188     {
1189       TopTools_IndexedDataMapOfShapeListOfShape solidsOfFace;
1190       if ( _sdVec[i]._hyp->IsToIgnoreShapes() )
1191         TopExp::MapShapesAndAncestors( _sdVec[i]._hypShape,
1192                                        TopAbs_FACE, TopAbs_SOLID, solidsOfFace);
1193
1194       exp.Init( _sdVec[i]._solid.Oriented( TopAbs_FORWARD ), TopAbs_FACE );
1195       for ( ; exp.More(); exp.Next() )
1196       {
1197         const TopoDS_Face& face = TopoDS::Face( exp.Current() );
1198         if ( helper.NbAncestors( face, *_mesh, TopAbs_SOLID ) < 2 )
1199           continue;
1200
1201         const TGeomID faceInd = getMeshDS()->ShapeToIndex( face );
1202         if ( _sdVec[i]._hyp->IsToIgnoreShapes() )
1203         {
1204           int nbSolids = solidsOfFace.FindFromKey( face ).Extent();
1205           if ( nbSolids > 1 )
1206             _sdVec[i]._ignoreFaceIds.insert( faceInd );
1207         }
1208
1209         if ( helper.IsReversedSubMesh( face ))
1210         {
1211           _sdVec[i]._reversedFaceIds.insert( faceInd );
1212         }
1213       }
1214     }
1215   }
1216
1217   // Find faces to shrink mesh on (solution 2 in issue 0020832);
1218   TopTools_IndexedMapOfShape shapes;
1219   for ( size_t i = 0; i < _sdVec.size(); ++i )
1220   {
1221     shapes.Clear();
1222     TopExp::MapShapes(_sdVec[i]._solid, TopAbs_EDGE, shapes);
1223     for ( int iE = 1; iE <= shapes.Extent(); ++iE )
1224     {
1225       const TopoDS_Shape& edge = shapes(iE);
1226       // find 2 faces sharing an edge
1227       TopoDS_Shape FF[2];
1228       PShapeIteratorPtr fIt = helper.GetAncestors(edge, *_mesh, TopAbs_FACE);
1229       while ( fIt->more())
1230       {
1231         const TopoDS_Shape* f = fIt->next();
1232         if ( helper.IsSubShape( *f, _sdVec[i]._solid))
1233           FF[ int( !FF[0].IsNull()) ] = *f;
1234       }
1235       if( FF[1].IsNull() ) continue; // seam edge can be shared by 1 FACE only
1236       // check presence of layers on them
1237       int ignore[2];
1238       for ( int j = 0; j < 2; ++j )
1239         ignore[j] = _sdVec[i]._ignoreFaceIds.count ( getMeshDS()->ShapeToIndex( FF[j] ));
1240       if ( ignore[0] == ignore[1] )
1241         continue; // nothing interesting
1242       TopoDS_Shape fWOL = FF[ ignore[0] ? 0 : 1 ];
1243       // check presence of layers on fWOL within an adjacent SOLID
1244       PShapeIteratorPtr sIt = helper.GetAncestors( fWOL, *_mesh, TopAbs_SOLID );
1245       while ( const TopoDS_Shape* solid = sIt->next() )
1246         if ( !solid->IsSame( _sdVec[i]._solid ))
1247         {
1248           int iSolid = solids.FindIndex( *solid );
1249           int  iFace = getMeshDS()->ShapeToIndex( fWOL );
1250           if ( iSolid > 0 && !_sdVec[ iSolid-1 ]._ignoreFaceIds.count( iFace ))
1251           {
1252             _sdVec[i]._noShrinkFaces.insert( iFace );
1253             fWOL.Nullify();
1254           }
1255         }
1256       // add edge to maps
1257       if ( !fWOL.IsNull())
1258       {
1259         TGeomID edgeInd = getMeshDS()->ShapeToIndex( edge );
1260         _sdVec[i]._shrinkShape2Shape.insert( make_pair( edgeInd, fWOL ));
1261       }
1262     }
1263   }
1264   // Exclude from _shrinkShape2Shape FACE's that can't be shrinked since
1265   // the algo of the SOLID sharing the FACE does not support it
1266   set< string > notSupportAlgos; notSupportAlgos.insert("Hexa_3D");
1267   for ( size_t i = 0; i < _sdVec.size(); ++i )
1268   {
1269     TopTools_MapOfShape noShrinkVertices;
1270     map< TGeomID, TopoDS_Shape >::iterator e2f = _sdVec[i]._shrinkShape2Shape.begin();
1271     for ( ; e2f != _sdVec[i]._shrinkShape2Shape.end(); ++e2f )
1272     {
1273       const TopoDS_Shape& fWOL = e2f->second;
1274       TGeomID           edgeID = e2f->first;
1275       bool notShrinkFace = false;
1276       PShapeIteratorPtr soIt = helper.GetAncestors(fWOL, *_mesh, TopAbs_SOLID);
1277       while ( soIt->more())
1278       {
1279         const TopoDS_Shape* solid = soIt->next();
1280         if ( _sdVec[i]._solid.IsSame( *solid )) continue;
1281         SMESH_Algo* algo = _mesh->GetGen()->GetAlgo( *_mesh, *solid );
1282         if ( !algo || !notSupportAlgos.count( algo->GetName() )) continue;
1283         notShrinkFace = true;
1284         for ( size_t j = 0; j < _sdVec.size(); ++j )
1285         {
1286           if ( _sdVec[j]._solid.IsSame( *solid ) )
1287             if ( _sdVec[j]._shrinkShape2Shape.count( edgeID ))
1288               notShrinkFace = false;
1289         }
1290       }
1291       if ( notShrinkFace )
1292       {
1293         _sdVec[i]._noShrinkFaces.insert( getMeshDS()->ShapeToIndex( fWOL ));
1294         for ( TopExp_Explorer vExp( fWOL, TopAbs_VERTEX ); vExp.More(); vExp.Next() )
1295           noShrinkVertices.Add( vExp.Current() );
1296       }
1297     }
1298     // erase from _shrinkShape2Shape all srink EDGE's of a SOLID connected
1299     // to the found not shrinked fWOL's
1300     e2f = _sdVec[i]._shrinkShape2Shape.begin();
1301     for ( ; e2f != _sdVec[i]._shrinkShape2Shape.end(); )
1302     {
1303       TGeomID edgeID = e2f->first;
1304       TopoDS_Vertex VV[2];
1305       TopExp::Vertices( TopoDS::Edge( getMeshDS()->IndexToShape( edgeID )),VV[0],VV[1]);
1306       if ( noShrinkVertices.Contains( VV[0] ) || noShrinkVertices.Contains( VV[1] ))
1307       {
1308         _sdVec[i]._noShrinkFaces.insert( getMeshDS()->ShapeToIndex( e2f->second ));
1309         _sdVec[i]._shrinkShape2Shape.erase( e2f++ );
1310       }
1311       else
1312       {
1313         e2f++;
1314       }
1315     }
1316   }
1317
1318   // Find the SHAPE along which to inflate _LayerEdge based on VERTEX
1319
1320   for ( size_t i = 0; i < _sdVec.size(); ++i )
1321   {
1322     shapes.Clear();
1323     TopExp::MapShapes(_sdVec[i]._solid, TopAbs_VERTEX, shapes);
1324     for ( int iV = 1; iV <= shapes.Extent(); ++iV )
1325     {
1326       const TopoDS_Shape& vertex = shapes(iV);
1327       // find faces WOL sharing the vertex
1328       vector< TopoDS_Shape > facesWOL;
1329       int totalNbFaces = 0;
1330       PShapeIteratorPtr fIt = helper.GetAncestors(vertex, *_mesh, TopAbs_FACE);
1331       while ( fIt->more())
1332       {
1333         const TopoDS_Shape* f = fIt->next();
1334         if ( helper.IsSubShape( *f, _sdVec[i]._solid ) )
1335         {
1336           totalNbFaces++;
1337           const int fID = getMeshDS()->ShapeToIndex( *f );
1338           if ( _sdVec[i]._ignoreFaceIds.count ( fID ) &&
1339                !_sdVec[i]._noShrinkFaces.count( fID ))
1340             facesWOL.push_back( *f );
1341         }
1342       }
1343       if ( facesWOL.size() == totalNbFaces || facesWOL.empty() )
1344         continue; // no layers at this vertex or no WOL
1345       TGeomID vInd = getMeshDS()->ShapeToIndex( vertex );
1346       switch ( facesWOL.size() )
1347       {
1348       case 1:
1349       {
1350         helper.SetSubShape( facesWOL[0] );
1351         if ( helper.IsRealSeam( vInd )) // inflate along a seam edge?
1352         {
1353           TopoDS_Shape seamEdge;
1354           PShapeIteratorPtr eIt = helper.GetAncestors(vertex, *_mesh, TopAbs_EDGE);
1355           while ( eIt->more() && seamEdge.IsNull() )
1356           {
1357             const TopoDS_Shape* e = eIt->next();
1358             if ( helper.IsRealSeam( *e ) )
1359               seamEdge = *e;
1360           }
1361           if ( !seamEdge.IsNull() )
1362           {
1363             _sdVec[i]._shrinkShape2Shape.insert( make_pair( vInd, seamEdge ));
1364             break;
1365           }
1366         }
1367         _sdVec[i]._shrinkShape2Shape.insert( make_pair( vInd, facesWOL[0] ));
1368         break;
1369       }
1370       case 2:
1371       {
1372         // find an edge shared by 2 faces
1373         PShapeIteratorPtr eIt = helper.GetAncestors(vertex, *_mesh, TopAbs_EDGE);
1374         while ( eIt->more())
1375         {
1376           const TopoDS_Shape* e = eIt->next();
1377           if ( helper.IsSubShape( *e, facesWOL[0]) &&
1378                helper.IsSubShape( *e, facesWOL[1]))
1379           {
1380             _sdVec[i]._shrinkShape2Shape.insert( make_pair( vInd, *e )); break;
1381           }
1382         }
1383         break;
1384       }
1385       default:
1386         return error("Not yet supported case", _sdVec[i]._index);
1387       }
1388     }
1389   }
1390
1391   return true;
1392 }
1393
1394 //================================================================================
1395 /*!
1396  * \brief Create the inner surface of the viscous layer and prepare data for infation
1397  */
1398 //================================================================================
1399
1400 bool _ViscousBuilder::makeLayer(_SolidData& data)
1401 {
1402   // get all sub-shapes to make layers on
1403   set<TGeomID> subIds, faceIds;
1404   subIds = data._noShrinkFaces;
1405   TopExp_Explorer exp( data._solid, TopAbs_FACE );
1406   for ( ; exp.More(); exp.Next() )
1407     {
1408       SMESH_subMesh* fSubM = _mesh->GetSubMesh( exp.Current() );
1409       if ( ! data._ignoreFaceIds.count( getMeshDS()->ShapeToIndex( exp.Current() )))
1410         faceIds.insert( fSubM->GetId() );
1411       SMESH_subMeshIteratorPtr subIt =
1412         fSubM->getDependsOnIterator(/*includeSelf=*/true, /*complexShapeFirst=*/false);
1413       while ( subIt->more() )
1414         subIds.insert( subIt->next()->GetId() );
1415     }
1416
1417   // make a map to find new nodes on sub-shapes shared with other SOLID
1418   map< TGeomID, TNode2Edge* > s2neMap;
1419   map< TGeomID, TNode2Edge* >::iterator s2ne;
1420   map< TGeomID, TopoDS_Shape >::iterator s2s = data._shrinkShape2Shape.begin();
1421   for (; s2s != data._shrinkShape2Shape.end(); ++s2s )
1422   {
1423     TGeomID shapeInd = s2s->first;
1424     for ( size_t i = 0; i < _sdVec.size(); ++i )
1425     {
1426       if ( _sdVec[i]._index == data._index ) continue;
1427       map< TGeomID, TopoDS_Shape >::iterator s2s2 = _sdVec[i]._shrinkShape2Shape.find( shapeInd );
1428       if ( s2s2 != _sdVec[i]._shrinkShape2Shape.end() &&
1429            *s2s == *s2s2 && !_sdVec[i]._n2eMap.empty() )
1430       {
1431         s2neMap.insert( make_pair( shapeInd, &_sdVec[i]._n2eMap ));
1432         break;
1433       }
1434     }
1435   }
1436
1437   // Create temporary faces and _LayerEdge's
1438
1439   dumpFunction(SMESH_Comment("makeLayers_")<<data._index); 
1440
1441   data._stepSize = Precision::Infinite();
1442   data._stepSizeNodes[0] = 0;
1443
1444   SMESH_MesherHelper helper( *_mesh );
1445   helper.SetSubShape( data._solid );
1446   helper.SetElementsOnShape(true);
1447
1448   vector< const SMDS_MeshNode*> newNodes; // of a mesh face
1449   TNode2Edge::iterator n2e2;
1450
1451   // collect _LayerEdge's of shapes they are based on
1452   const int nbShapes = getMeshDS()->MaxShapeIndex();
1453   vector< vector<_LayerEdge*> > edgesByGeom( nbShapes+1 );
1454
1455   for ( set<TGeomID>::iterator id = faceIds.begin(); id != faceIds.end(); ++id )
1456   {
1457     SMESHDS_SubMesh* smDS = getMeshDS()->MeshElements( *id );
1458     if ( !smDS ) return error(SMESH_Comment("Not meshed face ") << *id, data._index );
1459
1460     const TopoDS_Face& F = TopoDS::Face( getMeshDS()->IndexToShape( *id ));
1461     SMESH_ProxyMesh::SubMesh* proxySub =
1462       data._proxyMesh->getFaceSubM( F, /*create=*/true);
1463
1464     SMDS_ElemIteratorPtr eIt = smDS->GetElements();
1465     while ( eIt->more() )
1466     {
1467       const SMDS_MeshElement* face = eIt->next();
1468       newNodes.resize( face->NbCornerNodes() );
1469       double faceMaxCosin = -1;
1470       for ( int i = 0 ; i < face->NbCornerNodes(); ++i )
1471       {
1472         const SMDS_MeshNode* n = face->GetNode(i);
1473         TNode2Edge::iterator n2e = data._n2eMap.insert( make_pair( n, (_LayerEdge*)0 )).first;
1474         if ( !(*n2e).second )
1475         {
1476           // add a _LayerEdge
1477           _LayerEdge* edge = new _LayerEdge();
1478           n2e->second = edge;
1479           edge->_nodes.push_back( n );
1480           const int shapeID = n->getshapeId();
1481           edgesByGeom[ shapeID ].push_back( edge );
1482
1483           // set edge data or find already refined _LayerEdge and get data from it
1484           if ( n->GetPosition()->GetTypeOfPosition() != SMDS_TOP_FACE &&
1485                ( s2ne = s2neMap.find( shapeID )) != s2neMap.end() &&
1486                ( n2e2 = (*s2ne).second->find( n )) != s2ne->second->end())
1487           {
1488             _LayerEdge* foundEdge = (*n2e2).second;
1489             edge->Copy( *foundEdge, helper );
1490             // location of the last node is modified but we can restore
1491             // it by node position on _sWOL stored by the node
1492             const_cast< SMDS_MeshNode* >
1493               ( edge->_nodes.back() )->setXYZ( n->X(), n->Y(), n->Z() );
1494           }
1495           else
1496           {
1497             edge->_nodes.push_back( helper.AddNode( n->X(), n->Y(), n->Z() ));
1498             if ( !setEdgeData( *edge, subIds, helper, data ))
1499               return false;
1500           }
1501           dumpMove(edge->_nodes.back());
1502           if ( edge->_cosin > 0.01 )
1503           {
1504             if ( edge->_cosin > faceMaxCosin )
1505               faceMaxCosin = edge->_cosin;
1506           }
1507         }
1508         newNodes[ i ] = n2e->second->_nodes.back();
1509       }
1510       // create a temporary face
1511       const SMDS_MeshElement* newFace = new TmpMeshFace( newNodes, --_tmpFaceID );
1512       proxySub->AddElement( newFace );
1513
1514       // compute inflation step size by min size of element on a convex surface
1515       if ( faceMaxCosin > 0.1 )
1516         limitStepSize( data, face, faceMaxCosin );
1517     } // loop on 2D elements on a FACE
1518   } // loop on FACEs of a SOLID
1519
1520   data._epsilon = 1e-7;
1521   if ( data._stepSize < 1. )
1522     data._epsilon *= data._stepSize;
1523
1524   // Put _LayerEdge's into the vector data._edges
1525
1526   if ( !sortEdges( data, edgesByGeom ))
1527     return false;
1528
1529   // Set target nodes into _Simplex and _2NearEdges
1530   TNode2Edge::iterator n2e;
1531   for ( size_t i = 0; i < data._edges.size(); ++i )
1532   {
1533     if ( data._edges[i]->IsOnEdge())
1534       for ( int j = 0; j < 2; ++j )
1535       {
1536         if ( data._edges[i]->_nodes.back()->NbInverseElements(SMDSAbs_Volume) > 0 )
1537           break; // _LayerEdge is shared by two _SolidData's
1538         const SMDS_MeshNode* & n = data._edges[i]->_2neibors->_nodes[j];
1539         if (( n2e = data._n2eMap.find( n )) == data._n2eMap.end() )
1540           return error("_LayerEdge not found by src node", data._index);
1541         n = (*n2e).second->_nodes.back();
1542         data._edges[i]->_2neibors->_edges[j] = n2e->second;
1543       }
1544     else
1545       for ( size_t j = 0; j < data._edges[i]->_simplices.size(); ++j )
1546       {
1547         _Simplex& s = data._edges[i]->_simplices[j];
1548         s._nNext = data._n2eMap[ s._nNext ]->_nodes.back();
1549         s._nPrev = data._n2eMap[ s._nPrev ]->_nodes.back();
1550       }
1551   }
1552
1553   dumpFunctionEnd();
1554   return true;
1555 }
1556
1557 //================================================================================
1558 /*!
1559  * \brief Compute inflation step size by min size of element on a convex surface
1560  */
1561 //================================================================================
1562
1563 void _ViscousBuilder::limitStepSize( _SolidData&             data,
1564                                      const SMDS_MeshElement* face,
1565                                      const double            cosin)
1566 {
1567   int iN = 0;
1568   double minSize = 10 * data._stepSize;
1569   const int nbNodes = face->NbCornerNodes();
1570   for ( int i = 0; i < nbNodes; ++i )
1571   {
1572     const SMDS_MeshNode* nextN = face->GetNode( SMESH_MesherHelper::WrapIndex( i+1, nbNodes ));
1573     const SMDS_MeshNode* curN = face->GetNode( i );
1574     if ( nextN->GetPosition()->GetTypeOfPosition() == SMDS_TOP_FACE ||
1575          curN->GetPosition()->GetTypeOfPosition() == SMDS_TOP_FACE )
1576     {
1577       double dist = SMESH_TNodeXYZ( face->GetNode(i)).Distance( nextN );
1578       if ( dist < minSize )
1579         minSize = dist, iN = i;
1580     }
1581   }
1582   double newStep = 0.8 * minSize / cosin;
1583   if ( newStep < data._stepSize )
1584   {
1585     data._stepSize = newStep;
1586     data._stepSizeCoeff = 0.8 / cosin;
1587     data._stepSizeNodes[0] = face->GetNode( iN );
1588     data._stepSizeNodes[1] = face->GetNode( SMESH_MesherHelper::WrapIndex( iN+1, nbNodes ));
1589   }
1590 }
1591
1592 //================================================================================
1593 /*!
1594  * \brief Compute inflation step size by min size of element on a convex surface
1595  */
1596 //================================================================================
1597
1598 void _ViscousBuilder::limitStepSize( _SolidData& data, const double minSize)
1599 {
1600   if ( minSize < data._stepSize )
1601   {
1602     data._stepSize = minSize;
1603     if ( data._stepSizeNodes[0] )
1604     {
1605       double dist =
1606         SMESH_TNodeXYZ(data._stepSizeNodes[0]).Distance(data._stepSizeNodes[1]);
1607       data._stepSizeCoeff = data._stepSize / dist;
1608     }
1609   }
1610 }
1611
1612 //================================================================================
1613 /*!
1614  * \brief Separate shapes (and _LayerEdge's on them) to smooth from the rest ones
1615  */
1616 //================================================================================
1617
1618 bool _ViscousBuilder::sortEdges( _SolidData&                    data,
1619                                  vector< vector<_LayerEdge*> >& edgesByGeom)
1620 {
1621   // Find shapes needing smoothing; such a shape has _LayerEdge._normal on it's
1622   // boundry inclined at a sharp angle to the shape
1623
1624   list< TGeomID > shapesToSmooth;
1625   
1626   SMESH_MesherHelper helper( *_mesh );
1627   bool ok = true;
1628
1629   for ( size_t iS = 0; iS < edgesByGeom.size(); ++iS )
1630   {
1631     vector<_LayerEdge*>& eS = edgesByGeom[iS];
1632     if ( eS.empty() ) continue;
1633     TopoDS_Shape S = getMeshDS()->IndexToShape( iS );
1634     bool needSmooth = false;
1635     switch ( S.ShapeType() )
1636     {
1637     case TopAbs_EDGE: {
1638
1639       bool isShrinkEdge = !eS[0]->_sWOL.IsNull();
1640       for ( TopoDS_Iterator vIt( S ); vIt.More() && !needSmooth; vIt.Next() )
1641       {
1642         TGeomID iV = getMeshDS()->ShapeToIndex( vIt.Value() );
1643         vector<_LayerEdge*>& eV = edgesByGeom[ iV ];
1644         if ( eV.empty() ) continue;
1645         double cosin = eV[0]->_cosin;
1646         bool badCosin =
1647           ( !eV[0]->_sWOL.IsNull() && ( eV[0]->_sWOL.ShapeType() == TopAbs_EDGE || !isShrinkEdge));
1648         if ( badCosin )
1649         {
1650           gp_Vec dir1, dir2;
1651           if ( eV[0]->_sWOL.ShapeType() == TopAbs_EDGE )
1652             dir1 = getEdgeDir( TopoDS::Edge( eV[0]->_sWOL ), TopoDS::Vertex( vIt.Value() ));
1653           else
1654             dir1 = getFaceDir( TopoDS::Face( eV[0]->_sWOL ), TopoDS::Vertex( vIt.Value() ),
1655                                eV[0]->_nodes[0], helper, ok);
1656           dir2 = getEdgeDir( TopoDS::Edge( S ), TopoDS::Vertex( vIt.Value() ));
1657           double angle = dir1.Angle( dir2 );
1658           cosin = cos( angle );
1659         }
1660         needSmooth = ( cosin > 0.1 );
1661       }
1662       break;
1663     }
1664     case TopAbs_FACE: {
1665
1666       for ( TopExp_Explorer eExp( S, TopAbs_EDGE ); eExp.More() && !needSmooth; eExp.Next() )
1667       {
1668         TGeomID iE = getMeshDS()->ShapeToIndex( eExp.Current() );
1669         vector<_LayerEdge*>& eE = edgesByGeom[ iE ];
1670         if ( eE.empty() ) continue;
1671         if ( eE[0]->_sWOL.IsNull() )
1672         {
1673           for ( size_t i = 0; i < eE.size() && !needSmooth; ++i )
1674             needSmooth = ( eE[i]->_cosin > 0.1 );
1675         }
1676         else
1677         {
1678           const TopoDS_Face& F1 = TopoDS::Face( S );
1679           const TopoDS_Face& F2 = TopoDS::Face( eE[0]->_sWOL );
1680           const TopoDS_Edge& E  = TopoDS::Edge( eExp.Current() );
1681           for ( size_t i = 0; i < eE.size() && !needSmooth; ++i )
1682           {
1683             gp_Vec dir1 = getFaceDir( F1, E, eE[i]->_nodes[0], helper, ok );
1684             gp_Vec dir2 = getFaceDir( F2, E, eE[i]->_nodes[0], helper, ok );
1685             double angle = dir1.Angle( dir2 );
1686             double cosin = cos( angle );
1687             needSmooth = ( cosin > 0.1 );
1688           }
1689         }
1690       }
1691       break;
1692     }
1693     case TopAbs_VERTEX:
1694       continue;
1695     default:;
1696     }
1697     if ( needSmooth )
1698     {
1699       if ( S.ShapeType() == TopAbs_EDGE ) shapesToSmooth.push_front( iS );
1700       else                                shapesToSmooth.push_back ( iS );
1701     }
1702
1703   } // loop on edgesByGeom
1704
1705   data._edges.reserve( data._n2eMap.size() );
1706   data._endEdgeToSmooth.clear();
1707
1708   // first we put _LayerEdge's on shapes to smooth
1709   list< TGeomID >::iterator gIt = shapesToSmooth.begin();
1710   for ( ; gIt != shapesToSmooth.end(); ++gIt )
1711   {
1712     vector<_LayerEdge*>& eVec = edgesByGeom[ *gIt ];
1713     if ( eVec.empty() ) continue;
1714     data._edges.insert( data._edges.end(), eVec.begin(), eVec.end() );
1715     data._endEdgeToSmooth.push_back( data._edges.size() );
1716     eVec.clear();
1717   }
1718
1719   // then the rest _LayerEdge's
1720   for ( size_t iS = 0; iS < edgesByGeom.size(); ++iS )
1721   {
1722     vector<_LayerEdge*>& eVec = edgesByGeom[iS];
1723     data._edges.insert( data._edges.end(), eVec.begin(), eVec.end() );
1724     eVec.clear();
1725   }
1726
1727   return ok;
1728 }
1729
1730 //================================================================================
1731 /*!
1732  * \brief Set data of _LayerEdge needed for smoothing
1733  *  \param subIds - ids of sub-shapes of a SOLID to take into account faces from
1734  */
1735 //================================================================================
1736
1737 bool _ViscousBuilder::setEdgeData(_LayerEdge&         edge,
1738                                   const set<TGeomID>& subIds,
1739                                   SMESH_MesherHelper& helper,
1740                                   _SolidData&         data)
1741 {
1742   SMESH_MeshEditor editor(_mesh);
1743
1744   const SMDS_MeshNode* node = edge._nodes[0]; // source node
1745   SMDS_TypeOfPosition posType = node->GetPosition()->GetTypeOfPosition();
1746
1747   edge._len = 0;
1748   edge._2neibors = 0;
1749   edge._curvature = 0;
1750
1751   // --------------------------
1752   // Compute _normal and _cosin
1753   // --------------------------
1754
1755   edge._cosin = 0;
1756   edge._normal.SetCoord(0,0,0);
1757
1758   int totalNbFaces = 0;
1759   gp_Pnt p;
1760   gp_Vec du, dv, geomNorm;
1761   bool normOK = true;
1762
1763   TGeomID shapeInd = node->getshapeId();
1764   map< TGeomID, TopoDS_Shape >::const_iterator s2s = data._shrinkShape2Shape.find( shapeInd );
1765   bool onShrinkShape ( s2s != data._shrinkShape2Shape.end() );
1766   TopoDS_Shape vertEdge;
1767
1768   if ( onShrinkShape ) // one of faces the node is on has no layers
1769   {
1770     vertEdge = getMeshDS()->IndexToShape( s2s->first ); // vertex or edge
1771     if ( s2s->second.ShapeType() == TopAbs_EDGE )
1772     {
1773       // inflate from VERTEX along EDGE
1774       edge._normal = getEdgeDir( TopoDS::Edge( s2s->second ), TopoDS::Vertex( vertEdge ));
1775     }
1776     else if ( vertEdge.ShapeType() == TopAbs_VERTEX )
1777     {
1778       // inflate from VERTEX along FACE
1779       edge._normal = getFaceDir( TopoDS::Face( s2s->second ), TopoDS::Vertex( vertEdge ),
1780                                  node, helper, normOK, &edge._cosin);
1781     }
1782     else
1783     {
1784       // inflate from EDGE along FACE
1785       edge._normal = getFaceDir( TopoDS::Face( s2s->second ), TopoDS::Edge( vertEdge ),
1786                                  node, helper, normOK);
1787     }
1788   }
1789   else // layers are on all faces of SOLID the node is on
1790   {
1791     // find indices of geom faces the node lies on
1792     set<TGeomID> faceIds;
1793     if  ( posType == SMDS_TOP_FACE )
1794     {
1795       faceIds.insert( node->getshapeId() );
1796     }
1797     else
1798     {
1799       SMDS_ElemIteratorPtr fIt = node->GetInverseElementIterator(SMDSAbs_Face);
1800       while ( fIt->more() )
1801         faceIds.insert( editor.FindShape(fIt->next()));
1802     }
1803
1804     set<TGeomID>::iterator id = faceIds.begin();
1805     TopoDS_Face F;
1806     for ( ; id != faceIds.end(); ++id )
1807     {
1808       const TopoDS_Shape& s = getMeshDS()->IndexToShape( *id );
1809       if ( s.IsNull() || s.ShapeType() != TopAbs_FACE || !subIds.count( *id ))
1810         continue;
1811       totalNbFaces++;
1812       F = TopoDS::Face( s );
1813
1814       // IDEA: if there is a problem with finding a normal,
1815       // we can compute an area-weighted sum of normals of all faces sharing the node
1816       gp_XY uv = helper.GetNodeUV( F, node, 0, &normOK );
1817       Handle(Geom_Surface) surface = BRep_Tool::Surface( F );
1818       surface->D1( uv.X(), uv.Y(), p, du,dv );
1819       geomNorm = du ^ dv;
1820       double size2 = geomNorm.SquareMagnitude();
1821       if ( size2 < 1e-10 ) // singularity
1822       {
1823         SMDS_ElemIteratorPtr fIt = node->GetInverseElementIterator(SMDSAbs_Face);
1824         while ( fIt->more() )
1825         {
1826           const SMDS_MeshElement* f = fIt->next();
1827           if ( editor.FindShape( f ) == *id )
1828           {
1829             SMESH_MeshAlgos::FaceNormal( f, (gp_XYZ&) geomNorm.XYZ(), /*normalized=*/false );
1830             size2 = geomNorm.SquareMagnitude();
1831             break;
1832           }
1833         }
1834         // double ddu = 0, ddv = 0;
1835         // if ( du.SquareMagnitude() > dv.SquareMagnitude() )
1836         //   ddu = 1e-3;
1837         // else
1838         //   ddv = 1e-3;
1839         // surface->D1( uv.X()+ddu, uv.Y()+ddv, p, du,dv );
1840         // geomNorm = du ^ dv;
1841         // size2 = geomNorm.SquareMagnitude();
1842         // if ( size2 < 1e-10 )
1843         // {
1844         //   surface->D1( uv.X()-ddu, uv.Y()-ddv, p, du,dv );
1845         //   geomNorm = du ^ dv;
1846         //   size2 = geomNorm.SquareMagnitude();
1847         // }
1848       }
1849       if ( size2 > numeric_limits<double>::min() )
1850         geomNorm /= sqrt( size2 );
1851       else
1852         normOK = false;
1853       if ( helper.GetSubShapeOri( data._solid, F ) != TopAbs_REVERSED )
1854         geomNorm.Reverse();
1855       edge._normal += geomNorm.XYZ();
1856     }
1857     if ( totalNbFaces == 0 )
1858       return error(SMESH_Comment("Can't get normal to node ") << node->GetID(), data._index);
1859
1860     edge._normal /= totalNbFaces;
1861
1862     switch ( posType )
1863     {
1864     case SMDS_TOP_FACE:
1865       edge._cosin = 0; break;
1866
1867     case SMDS_TOP_EDGE: {
1868       TopoDS_Edge E = TopoDS::Edge( helper.GetSubShapeByNode( node, getMeshDS()));
1869       gp_Vec inFaceDir = getFaceDir( F, E, node, helper, normOK);
1870       double angle = inFaceDir.Angle( edge._normal ); // [0,PI]
1871       edge._cosin = cos( angle );
1872       //cout << "Cosin on EDGE " << edge._cosin << " node " << node->GetID() << endl;
1873       break;
1874     }
1875     case SMDS_TOP_VERTEX: {
1876       TopoDS_Vertex V = TopoDS::Vertex( helper.GetSubShapeByNode( node, getMeshDS()));
1877       gp_XYZ inFaceDir = getFaceDir( F, V, node, helper, normOK);
1878       double angle = gp_Vec( inFaceDir).Angle( edge._normal ); // [0,PI]
1879       edge._cosin = cos( angle );
1880       //cout << "Cosin on VERTEX " << edge._cosin << " node " << node->GetID() << endl;
1881       break;
1882     }
1883     default:
1884       return error(SMESH_Comment("Invalid shape position of node ")<<node, data._index);
1885     }
1886   }
1887
1888   double normSize = edge._normal.SquareModulus();
1889   if ( normSize < numeric_limits<double>::min() )
1890     return error(SMESH_Comment("Bad normal at node ")<< node->GetID(), data._index );
1891
1892   edge._normal /= sqrt( normSize );
1893
1894   // TODO: if ( !normOK ) then get normal by mesh faces
1895
1896   // Set the rest data
1897   // --------------------
1898   if ( onShrinkShape )
1899   {
1900     edge._sWOL = (*s2s).second;
1901
1902     SMDS_MeshNode* tgtNode = const_cast<SMDS_MeshNode*>( edge._nodes.back() );
1903     if ( SMESHDS_SubMesh* sm = getMeshDS()->MeshElements( data._solid ))
1904       sm->RemoveNode( tgtNode , /*isNodeDeleted=*/false );
1905
1906     // set initial position which is parameters on _sWOL in this case
1907     if ( edge._sWOL.ShapeType() == TopAbs_EDGE )
1908     {
1909       double u = helper.GetNodeU( TopoDS::Edge( edge._sWOL ), node, 0, &normOK );
1910       edge._pos.push_back( gp_XYZ( u, 0, 0));
1911       getMeshDS()->SetNodeOnEdge( tgtNode, TopoDS::Edge( edge._sWOL ), u );
1912     }
1913     else // TopAbs_FACE
1914     {
1915       gp_XY uv = helper.GetNodeUV( TopoDS::Face( edge._sWOL ), node, 0, &normOK );
1916       edge._pos.push_back( gp_XYZ( uv.X(), uv.Y(), 0));
1917       getMeshDS()->SetNodeOnFace( tgtNode, TopoDS::Face( edge._sWOL ), uv.X(), uv.Y() );
1918     }
1919   }
1920   else
1921   {
1922     edge._pos.push_back( SMESH_TNodeXYZ( node ));
1923
1924     if ( posType == SMDS_TOP_FACE )
1925     {
1926       getSimplices( node, edge._simplices, data._ignoreFaceIds, &data );
1927       double avgNormProj = 0, avgLen = 0;
1928       for ( size_t i = 0; i < edge._simplices.size(); ++i )
1929       {
1930         gp_XYZ vec = edge._pos.back() - SMESH_TNodeXYZ( edge._simplices[i]._nPrev );
1931         avgNormProj += edge._normal * vec;
1932         avgLen += vec.Modulus();
1933       }
1934       avgNormProj /= edge._simplices.size();
1935       avgLen /= edge._simplices.size();
1936       edge._curvature = _Curvature::New( avgNormProj, avgLen );
1937     }
1938   }
1939
1940   // Set neighbour nodes for a _LayerEdge based on EDGE
1941
1942   if ( posType == SMDS_TOP_EDGE /*||
1943        ( onShrinkShape && posType == SMDS_TOP_VERTEX && fabs( edge._cosin ) < 1e-10 )*/)
1944   {
1945     edge._2neibors = new _2NearEdges;
1946     // target node instead of source ones will be set later
1947     if ( ! findNeiborsOnEdge( &edge,
1948                               edge._2neibors->_nodes[0],
1949                               edge._2neibors->_nodes[1],
1950                               data))
1951       return false;
1952     edge.SetDataByNeighbors( edge._2neibors->_nodes[0],
1953                              edge._2neibors->_nodes[1],
1954                              helper);
1955   }
1956
1957   edge.SetCosin( edge._cosin ); // to update edge._lenFactor
1958
1959   return true;
1960 }
1961
1962 //================================================================================
1963 /*!
1964  * \brief Find 2 neigbor nodes of a node on EDGE
1965  */
1966 //================================================================================
1967
1968 bool _ViscousBuilder::findNeiborsOnEdge(const _LayerEdge*     edge,
1969                                         const SMDS_MeshNode*& n1,
1970                                         const SMDS_MeshNode*& n2,
1971                                         _SolidData&           data)
1972 {
1973   const SMDS_MeshNode* node = edge->_nodes[0];
1974   const int shapeInd = node->getshapeId();
1975   SMESHDS_SubMesh* edgeSM = 0;
1976   if ( node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_EDGE )
1977   {
1978     
1979     edgeSM = getMeshDS()->MeshElements( shapeInd );
1980     if ( !edgeSM || edgeSM->NbElements() == 0 )
1981       return error(SMESH_Comment("Not meshed EDGE ") << shapeInd, data._index);
1982   }
1983   int iN = 0;
1984   n2 = 0;
1985   SMDS_ElemIteratorPtr eIt = node->GetInverseElementIterator(SMDSAbs_Edge);
1986   while ( eIt->more() && !n2 )
1987   {
1988     const SMDS_MeshElement* e = eIt->next();
1989     const SMDS_MeshNode*   nNeibor = e->GetNode( 0 );
1990     if ( nNeibor == node ) nNeibor = e->GetNode( 1 );
1991     if ( edgeSM )
1992     {
1993       if (!edgeSM->Contains(e)) continue;
1994     }
1995     else
1996     {
1997       TopoDS_Shape s = SMESH_MesherHelper::GetSubShapeByNode(nNeibor, getMeshDS() );
1998       if ( !SMESH_MesherHelper::IsSubShape( s, edge->_sWOL )) continue;
1999     }
2000     ( iN++ ? n2 : n1 ) = nNeibor;
2001   }
2002   if ( !n2 )
2003     return error(SMESH_Comment("Wrongly meshed EDGE ") << shapeInd, data._index);
2004   return true;
2005 }
2006
2007 //================================================================================
2008 /*!
2009  * \brief Set _curvature and _2neibors->_plnNorm by 2 neigbor nodes residing the same EDGE
2010  */
2011 //================================================================================
2012
2013 void _LayerEdge::SetDataByNeighbors( const SMDS_MeshNode* n1,
2014                                      const SMDS_MeshNode* n2,
2015                                      SMESH_MesherHelper&  helper)
2016 {
2017   if ( _nodes[0]->GetPosition()->GetTypeOfPosition() != SMDS_TOP_EDGE )
2018     return;
2019
2020   gp_XYZ pos = SMESH_TNodeXYZ( _nodes[0] );
2021   gp_XYZ vec1 = pos - SMESH_TNodeXYZ( n1 );
2022   gp_XYZ vec2 = pos - SMESH_TNodeXYZ( n2 );
2023
2024   // Set _curvature
2025
2026   double sumLen = vec1.Modulus() + vec2.Modulus();
2027   _2neibors->_wgt[0] = 1 - vec1.Modulus() / sumLen;
2028   _2neibors->_wgt[1] = 1 - vec2.Modulus() / sumLen;
2029   double avgNormProj = 0.5 * ( _normal * vec1 + _normal * vec2 );
2030   double avgLen = 0.5 * ( vec1.Modulus() + vec2.Modulus() );
2031   if ( _curvature ) delete _curvature;
2032   _curvature = _Curvature::New( avgNormProj, avgLen );
2033 #ifdef __myDEBUG
2034 //     if ( _curvature )
2035 //       cout << _nodes[0]->GetID()
2036 //            << " CURV r,k: " << _curvature->_r<<","<<_curvature->_k
2037 //            << " proj = "<<avgNormProj<< " len = " << avgLen << "| lenDelta(0) = "
2038 //            << _curvature->lenDelta(0) << endl;
2039 #endif
2040
2041   // Set _plnNorm
2042
2043   if ( _sWOL.IsNull() )
2044   {
2045     TopoDS_Shape S = helper.GetSubShapeByNode( _nodes[0], helper.GetMeshDS() );
2046     gp_XYZ dirE = getEdgeDir( TopoDS::Edge( S ), _nodes[0], helper );
2047     gp_XYZ plnNorm = dirE ^ _normal;
2048     double proj0 = plnNorm * vec1;
2049     double proj1 = plnNorm * vec2;
2050     if ( fabs( proj0 ) > 1e-10 || fabs( proj1 ) > 1e-10 )
2051     {
2052       if ( _2neibors->_plnNorm ) delete _2neibors->_plnNorm;
2053       _2neibors->_plnNorm = new gp_XYZ( plnNorm.Normalized() );
2054     }
2055   }
2056 }
2057
2058 //================================================================================
2059 /*!
2060  * \brief Copy data from a _LayerEdge of other SOLID and based on the same node;
2061  * this and other _LayerEdge's are inflated along a FACE or an EDGE
2062  */
2063 //================================================================================
2064
2065 void _LayerEdge::Copy( _LayerEdge& other, SMESH_MesherHelper& helper )
2066 {
2067   _nodes     = other._nodes;
2068   _normal    = other._normal;
2069   _len       = 0;
2070   _lenFactor = other._lenFactor;
2071   _cosin     = other._cosin;
2072   _sWOL      = other._sWOL;
2073   _2neibors  = other._2neibors;
2074   _curvature = 0; std::swap( _curvature, other._curvature );
2075   _2neibors  = 0; std::swap( _2neibors,  other._2neibors );
2076
2077   if ( _sWOL.ShapeType() == TopAbs_EDGE )
2078   {
2079     double u = helper.GetNodeU( TopoDS::Edge( _sWOL ), _nodes[0] );
2080     _pos.push_back( gp_XYZ( u, 0, 0));
2081   }
2082   else // TopAbs_FACE
2083   {
2084     gp_XY uv = helper.GetNodeUV( TopoDS::Face( _sWOL ), _nodes[0]);
2085     _pos.push_back( gp_XYZ( uv.X(), uv.Y(), 0));
2086   }
2087 }
2088
2089 //================================================================================
2090 /*!
2091  * \brief Set _cosin and _lenFactor
2092  */
2093 //================================================================================
2094
2095 void _LayerEdge::SetCosin( double cosin )
2096 {
2097   _cosin = cosin;
2098   _lenFactor = ( _cosin > 0.1 ) ?  1./sqrt(1-_cosin*_cosin) : 1.0;
2099 }
2100
2101 //================================================================================
2102 /*!
2103  * \brief Fills a vector<_Simplex > 
2104  */
2105 //================================================================================
2106
2107 void _ViscousBuilder::getSimplices( const SMDS_MeshNode* node,
2108                                     vector<_Simplex>&    simplices,
2109                                     const set<TGeomID>&  ingnoreShapes,
2110                                     const _SolidData*    dataToCheckOri,
2111                                     const bool           toSort)
2112 {
2113   simplices.clear();
2114   SMDS_ElemIteratorPtr fIt = node->GetInverseElementIterator(SMDSAbs_Face);
2115   while ( fIt->more() )
2116   {
2117     const SMDS_MeshElement* f = fIt->next();
2118     const TGeomID    shapeInd = f->getshapeId();
2119     if ( ingnoreShapes.count( shapeInd )) continue;
2120     const int nbNodes = f->NbCornerNodes();
2121     const int  srcInd = f->GetNodeIndex( node );
2122     const SMDS_MeshNode* nPrev = f->GetNode( SMESH_MesherHelper::WrapIndex( srcInd-1, nbNodes ));
2123     const SMDS_MeshNode* nNext = f->GetNode( SMESH_MesherHelper::WrapIndex( srcInd+1, nbNodes ));
2124     const SMDS_MeshNode* nOpp  = f->GetNode( SMESH_MesherHelper::WrapIndex( srcInd+2, nbNodes ));
2125     if ( dataToCheckOri && dataToCheckOri->_reversedFaceIds.count( shapeInd ))
2126       std::swap( nPrev, nNext );
2127     simplices.push_back( _Simplex( nPrev, nNext, nOpp ));
2128   }
2129
2130   if ( toSort )
2131   {
2132     vector<_Simplex> sortedSimplices( simplices.size() );
2133     sortedSimplices[0] = simplices[0];
2134     int nbFound = 0;
2135     for ( size_t i = 1; i < simplices.size(); ++i )
2136     {
2137       for ( size_t j = 1; j < simplices.size(); ++j )
2138         if ( sortedSimplices[i-1]._nNext == simplices[j]._nPrev )
2139         {
2140           sortedSimplices[i] = simplices[j];
2141           nbFound++;
2142           break;
2143         }
2144     }
2145     if ( nbFound == simplices.size() - 1 )
2146       simplices.swap( sortedSimplices );
2147   }
2148 }
2149
2150 //================================================================================
2151 /*!
2152  * \brief DEBUG. Create groups contating temorary data of _LayerEdge's
2153  */
2154 //================================================================================
2155
2156 void _ViscousBuilder::makeGroupOfLE()
2157 {
2158 #ifdef _DEBUG_
2159   for ( size_t i = 0 ; i < _sdVec.size(); ++i )
2160   {
2161     if ( _sdVec[i]._edges.empty() ) continue;
2162 //     string name = SMESH_Comment("_LayerEdge's_") << i;
2163 //     int id;
2164 //     SMESH_Group* g = _mesh->AddGroup(SMDSAbs_Edge, name.c_str(), id );
2165 //     SMESHDS_Group* gDS = (SMESHDS_Group*)g->GetGroupDS();
2166 //     SMESHDS_Mesh* mDS = _mesh->GetMeshDS();
2167
2168     dumpFunction( SMESH_Comment("make_LayerEdge_") << i );
2169     for ( size_t j = 0 ; j < _sdVec[i]._edges.size(); ++j )
2170     {
2171       _LayerEdge* le = _sdVec[i]._edges[j];
2172       for ( size_t iN = 1; iN < le->_nodes.size(); ++iN )
2173         dumpCmd(SMESH_Comment("mesh.AddEdge([ ") <<le->_nodes[iN-1]->GetID()
2174                 << ", " << le->_nodes[iN]->GetID() <<"])");
2175       //gDS->SMDSGroup().Add( mDS->AddEdge( le->_nodes[iN-1], le->_nodes[iN]));
2176     }
2177     dumpFunctionEnd();
2178
2179     dumpFunction( SMESH_Comment("makeNormals") << i );
2180     for ( size_t j = 0 ; j < _sdVec[i]._edges.size(); ++j )
2181     {
2182       _LayerEdge& edge = *_sdVec[i]._edges[j];
2183       SMESH_TNodeXYZ nXYZ( edge._nodes[0] );
2184       nXYZ += edge._normal * _sdVec[i]._stepSize;
2185       dumpCmd(SMESH_Comment("mesh.AddEdge([ ") <<edge._nodes[0]->GetID()
2186               << ", mesh.AddNode( " << nXYZ.X()<<","<< nXYZ.Y()<<","<< nXYZ.Z()<<")])");
2187     }
2188     dumpFunctionEnd();
2189
2190 //     name = SMESH_Comment("tmp_faces ") << i;
2191 //     g = _mesh->AddGroup(SMDSAbs_Face, name.c_str(), id );
2192 //     gDS = (SMESHDS_Group*)g->GetGroupDS();
2193 //     SMESH_MeshEditor editor( _mesh );
2194     dumpFunction( SMESH_Comment("makeTmpFaces_") << i );
2195     TopExp_Explorer fExp( _sdVec[i]._solid, TopAbs_FACE );
2196     for ( ; fExp.More(); fExp.Next() )
2197     {
2198       if (const SMESHDS_SubMesh* sm = _sdVec[i]._proxyMesh->GetProxySubMesh( fExp.Current()))
2199       {
2200         SMDS_ElemIteratorPtr fIt = sm->GetElements();
2201         while ( fIt->more())
2202         {
2203           const SMDS_MeshElement* e = fIt->next();
2204           SMESH_Comment cmd("mesh.AddFace([");
2205           for ( int j=0; j < e->NbCornerNodes(); ++j )
2206             cmd << e->GetNode(j)->GetID() << (j+1<e->NbCornerNodes() ? ",": "])");
2207           dumpCmd( cmd );
2208           //vector<const SMDS_MeshNode*> nodes( e->begin_nodes(), e->end_nodes() );
2209           //gDS->SMDSGroup().Add( editor.AddElement( nodes, e->GetType(), e->IsPoly()));
2210         }
2211       }
2212     }
2213     dumpFunctionEnd();
2214   }
2215 #endif
2216 }
2217
2218 //================================================================================
2219 /*!
2220  * \brief Increase length of _LayerEdge's to reach the required thickness of layers
2221  */
2222 //================================================================================
2223
2224 bool _ViscousBuilder::inflate(_SolidData& data)
2225 {
2226   SMESH_MesherHelper helper( *_mesh );
2227
2228   // Limit inflation step size by geometry size found by itersecting
2229   // normals of _LayerEdge's with mesh faces
2230   double geomSize = Precision::Infinite(), intersecDist;
2231   auto_ptr<SMESH_ElementSearcher> searcher
2232     ( SMESH_MeshAlgos::GetElementSearcher( *getMeshDS(),
2233                                            data._proxyMesh->GetFaces( data._solid )) );
2234   for ( size_t i = 0; i < data._edges.size(); ++i )
2235   {
2236     if ( data._edges[i]->IsOnEdge() ) continue;
2237     data._edges[i]->FindIntersection( *searcher, intersecDist, data._epsilon );
2238     if ( geomSize > intersecDist && intersecDist > 0 )
2239       geomSize = intersecDist;
2240   }
2241   if ( data._stepSize > 0.3 * geomSize )
2242     limitStepSize( data, 0.3 * geomSize );
2243
2244   const double tgtThick = data._hyp->GetTotalThickness();
2245   if ( data._stepSize > tgtThick )
2246     limitStepSize( data, tgtThick );
2247
2248   if ( data._stepSize < 1. )
2249     data._epsilon = data._stepSize * 1e-7;
2250
2251 #ifdef __myDEBUG
2252   cout << "-- geomSize = " << geomSize << ", stepSize = " << data._stepSize << endl;
2253 #endif
2254
2255   double avgThick = 0, curThick = 0, distToIntersection = Precision::Infinite();
2256   int nbSteps = 0, nbRepeats = 0;
2257   while ( 1.01 * avgThick < tgtThick )
2258   {
2259     // new target length
2260     curThick += data._stepSize;
2261     if ( curThick > tgtThick )
2262     {
2263       curThick = tgtThick + ( tgtThick-avgThick ) * nbRepeats;
2264       nbRepeats++;
2265     }
2266
2267     // Elongate _LayerEdge's
2268     dumpFunction(SMESH_Comment("inflate")<<data._index<<"_step"<<nbSteps); // debug
2269     for ( size_t i = 0; i < data._edges.size(); ++i )
2270     {
2271       data._edges[i]->SetNewLength( curThick, helper );
2272     }
2273     dumpFunctionEnd();
2274
2275     if ( !nbSteps )
2276       if ( !updateNormals( data, helper ) )
2277         return false;
2278
2279     // Improve and check quality
2280     if ( !smoothAndCheck( data, nbSteps, distToIntersection ))
2281     {
2282       if ( nbSteps > 0 )
2283       {
2284         dumpFunction(SMESH_Comment("invalidate")<<data._index<<"_step"<<nbSteps); // debug
2285         for ( size_t i = 0; i < data._edges.size(); ++i )
2286         {
2287           data._edges[i]->InvalidateStep( nbSteps+1 );
2288         }
2289         dumpFunctionEnd();
2290       }
2291       break; // no more inflating possible
2292     }
2293     nbSteps++;
2294
2295     // Evaluate achieved thickness
2296     avgThick = 0;
2297     for ( size_t i = 0; i < data._edges.size(); ++i )
2298       avgThick += data._edges[i]->_len;
2299     avgThick /= data._edges.size();
2300 #ifdef __myDEBUG
2301     cout << "-- Thickness " << avgThick << " reached" << endl;
2302 #endif
2303
2304     if ( distToIntersection < avgThick*1.5 )
2305     {
2306 #ifdef __myDEBUG
2307       cout << "-- Stop inflation since distToIntersection( "<<distToIntersection<<" ) < avgThick( "
2308            << avgThick << " ) * 1.5" << endl;
2309 #endif
2310       break;
2311     }
2312     // new step size
2313     limitStepSize( data, 0.25 * distToIntersection );
2314     if ( data._stepSizeNodes[0] )
2315       data._stepSize = data._stepSizeCoeff *
2316         SMESH_TNodeXYZ(data._stepSizeNodes[0]).Distance(data._stepSizeNodes[1]);
2317   }
2318
2319   if (nbSteps == 0 )
2320     return error("failed at the very first inflation step", data._index);
2321
2322   return true;
2323 }
2324
2325 //================================================================================
2326 /*!
2327  * \brief Improve quality of layer inner surface and check intersection
2328  */
2329 //================================================================================
2330
2331 bool _ViscousBuilder::smoothAndCheck(_SolidData& data,
2332                                      const int   nbSteps,
2333                                      double &    distToIntersection)
2334 {
2335   if ( data._endEdgeToSmooth.empty() )
2336     return true; // no shapes needing smoothing
2337
2338   bool moved, improved;
2339
2340   SMESH_MesherHelper helper(*_mesh);
2341   Handle(Geom_Surface) surface;
2342   TopoDS_Face F;
2343
2344   int iBeg, iEnd = 0;
2345   for ( size_t iS = 0; iS < data._endEdgeToSmooth.size(); ++iS )
2346   {
2347     iBeg = iEnd;
2348     iEnd = data._endEdgeToSmooth[ iS ];
2349
2350     if ( !data._edges[ iBeg ]->_sWOL.IsNull() &&
2351          data._edges[ iBeg ]->_sWOL.ShapeType() == TopAbs_FACE )
2352     {
2353       if ( !F.IsSame( data._edges[ iBeg ]->_sWOL )) {
2354         F = TopoDS::Face( data._edges[ iBeg ]->_sWOL );
2355         helper.SetSubShape( F );
2356         surface = BRep_Tool::Surface( F );
2357       }
2358     }
2359     else
2360     {
2361       F.Nullify(); surface.Nullify();
2362     }
2363     TGeomID sInd = data._edges[ iBeg ]->_nodes[0]->getshapeId();
2364
2365     if ( data._edges[ iBeg ]->IsOnEdge() )
2366     { 
2367       dumpFunction(SMESH_Comment("smooth")<<data._index << "_Ed"<<sInd <<"_InfStep"<<nbSteps);
2368
2369       // try a simple solution on an analytic EDGE
2370       if ( !smoothAnalyticEdge( data, iBeg, iEnd, surface, F, helper ))
2371       {
2372         // smooth on EDGE's
2373         int step = 0;
2374         do {
2375           moved = false;
2376           for ( int i = iBeg; i < iEnd; ++i )
2377           {
2378             moved |= data._edges[i]->SmoothOnEdge(surface, F, helper);
2379           }
2380           dumpCmd( SMESH_Comment("# end step ")<<step);
2381         }
2382         while ( moved && step++ < 5 );
2383         //cout << " NB STEPS: " << step << endl;
2384       }
2385       dumpFunctionEnd();
2386     }
2387     else
2388     {
2389       // smooth on FACE's
2390       int step = 0, badNb = 0; moved = true;
2391       while (( ++step <= 5 && moved ) || improved )
2392       {
2393         dumpFunction(SMESH_Comment("smooth")<<data._index<<"_Fa"<<sInd
2394                      <<"_InfStep"<<nbSteps<<"_"<<step); // debug
2395         int oldBadNb = badNb;
2396         badNb = 0;
2397         moved = false;
2398         for ( int i = iBeg; i < iEnd; ++i )
2399           moved |= data._edges[i]->Smooth(badNb);
2400         improved = ( badNb < oldBadNb );
2401
2402         dumpFunctionEnd();
2403       }
2404       if ( badNb > 0 )
2405       {
2406 #ifdef __myDEBUG
2407         for ( int i = iBeg; i < iEnd; ++i )
2408         {
2409           _LayerEdge* edge = data._edges[i];
2410           SMESH_TNodeXYZ tgtXYZ( edge->_nodes.back() );
2411           for ( size_t j = 0; j < edge->_simplices.size(); ++j )
2412             if ( !edge->_simplices[j].IsForward( edge->_nodes[0], &tgtXYZ ))
2413             {
2414               cout << "Bad simplex ( " << edge->_nodes[0]->GetID()<< " "<< tgtXYZ._node->GetID()
2415                    << " "<< edge->_simplices[j]._nPrev->GetID()
2416                    << " "<< edge->_simplices[j]._nNext->GetID() << " )" << endl;
2417               return false;
2418             }
2419         }
2420 #endif
2421         return false;
2422       }
2423     }
2424   } // loop on shapes to smooth
2425
2426   // Check if the last segments of _LayerEdge intersects 2D elements;
2427   // checked elements are either temporary faces or faces on surfaces w/o the layers
2428
2429   auto_ptr<SMESH_ElementSearcher> searcher
2430     ( SMESH_MeshAlgos::GetElementSearcher( *getMeshDS(),
2431                                            data._proxyMesh->GetFaces( data._solid )) );
2432
2433   distToIntersection = Precision::Infinite();
2434   double dist;
2435   const SMDS_MeshElement* intFace = 0;
2436 #ifdef __myDEBUG
2437   const SMDS_MeshElement* closestFace = 0;
2438   int iLE = 0;
2439 #endif
2440   for ( size_t i = 0; i < data._edges.size(); ++i )
2441   {
2442     if ( data._edges[i]->FindIntersection( *searcher, dist, data._epsilon, &intFace ))
2443       return false;
2444     if ( distToIntersection > dist )
2445     {
2446       distToIntersection = dist;
2447 #ifdef __myDEBUG
2448       iLE = i;
2449       closestFace = intFace;
2450 #endif
2451     }
2452   }
2453 #ifdef __myDEBUG
2454   if ( closestFace )
2455   {
2456     SMDS_MeshElement::iterator nIt = closestFace->begin_nodes();
2457     cout << "Shortest distance: _LayerEdge nodes: tgt " << data._edges[iLE]->_nodes.back()->GetID()
2458          << " src " << data._edges[iLE]->_nodes[0]->GetID()<< ", intersection with face ("
2459          << (*nIt++)->GetID()<<" "<< (*nIt++)->GetID()<<" "<< (*nIt++)->GetID()
2460          << ") distance = " << distToIntersection<< endl;
2461   }
2462 #endif
2463
2464   return true;
2465 }
2466
2467 //================================================================================
2468 /*!
2469  * \brief Return a curve of the EDGE to be used for smoothing and arrange
2470  *        _LayerEdge's to be in a consequent order
2471  */
2472 //================================================================================
2473
2474 Handle(Geom_Curve) _SolidData::CurveForSmooth( const TopoDS_Edge&    E,
2475                                                const int             iFrom,
2476                                                const int             iTo,
2477                                                Handle(Geom_Surface)& surface,
2478                                                const TopoDS_Face&    F,
2479                                                SMESH_MesherHelper&   helper)
2480 {
2481   TGeomID eIndex = helper.GetMeshDS()->ShapeToIndex( E );
2482
2483   map< TGeomID, Handle(Geom_Curve)>::iterator i2curve = _edge2curve.find( eIndex );
2484
2485   if ( i2curve == _edge2curve.end() )
2486   {
2487     // sort _LayerEdge's by position on the EDGE
2488     {
2489       map< double, _LayerEdge* > u2edge;
2490       for ( int i = iFrom; i < iTo; ++i )
2491         u2edge.insert( make_pair( helper.GetNodeU( E, _edges[i]->_nodes[0] ), _edges[i] ));
2492
2493       ASSERT( u2edge.size() == iTo - iFrom );
2494       map< double, _LayerEdge* >::iterator u2e = u2edge.begin();
2495       for ( int i = iFrom; i < iTo; ++i, ++u2e )
2496         _edges[i] = u2e->second;
2497
2498       // set _2neibors according to the new order
2499       for ( int i = iFrom; i < iTo-1; ++i )
2500         if ( _edges[i]->_2neibors->_nodes[1] != _edges[i+1]->_nodes.back() )
2501           _edges[i]->_2neibors->reverse();
2502       if ( u2edge.size() > 1 &&
2503            _edges[iTo-1]->_2neibors->_nodes[0] != _edges[iTo-2]->_nodes.back() )
2504         _edges[iTo-1]->_2neibors->reverse();
2505     }
2506
2507     SMESHDS_SubMesh* smDS = helper.GetMeshDS()->MeshElements( eIndex );
2508
2509     TopLoc_Location loc; double f,l;
2510
2511     Handle(Geom_Line)   line;
2512     Handle(Geom_Circle) circle;
2513     bool isLine, isCirc;
2514     if ( F.IsNull() ) // 3D case
2515     {
2516       // check if the EDGE is a line
2517       Handle(Geom_Curve) curve = BRep_Tool::Curve( E, loc, f, l);
2518       if ( curve->IsKind( STANDARD_TYPE( Geom_TrimmedCurve )))
2519         curve = Handle(Geom_TrimmedCurve)::DownCast( curve )->BasisCurve();
2520
2521       line   = Handle(Geom_Line)::DownCast( curve );
2522       circle = Handle(Geom_Circle)::DownCast( curve );
2523       isLine = (!line.IsNull());
2524       isCirc = (!circle.IsNull());
2525
2526       if ( !isLine && !isCirc ) // Check if the EDGE is close to a line
2527       {
2528         Bnd_B3d bndBox;
2529         SMDS_NodeIteratorPtr nIt = smDS->GetNodes();
2530         while ( nIt->more() )
2531           bndBox.Add( SMESH_TNodeXYZ( nIt->next() ));
2532         gp_XYZ size = bndBox.CornerMax() - bndBox.CornerMin();
2533
2534         SMESH_TNodeXYZ p0( _edges[iFrom]->_2neibors->_nodes[0] );
2535         SMESH_TNodeXYZ p1( _edges[iFrom]->_2neibors->_nodes[1] );
2536         const double lineTol = 1e-2 * ( p0 - p1 ).Modulus();
2537         for ( int i = 0; i < 3 && !isLine; ++i )
2538           isLine = ( size.Coord( i+1 ) <= lineTol );
2539       }
2540       if ( !isLine && !isCirc && iTo-iFrom > 2) // Check if the EDGE is close to a circle
2541       {
2542         // TODO
2543       }
2544     }
2545     else // 2D case
2546     {
2547       // check if the EDGE is a line
2548       Handle(Geom2d_Curve) curve = BRep_Tool::CurveOnSurface( E, F, f, l);
2549       if ( curve->IsKind( STANDARD_TYPE( Geom2d_TrimmedCurve )))
2550         curve = Handle(Geom2d_TrimmedCurve)::DownCast( curve )->BasisCurve();
2551
2552       Handle(Geom2d_Line)   line2d   = Handle(Geom2d_Line)::DownCast( curve );
2553       Handle(Geom2d_Circle) circle2d = Handle(Geom2d_Circle)::DownCast( curve );
2554       isLine = (!line2d.IsNull());
2555       isCirc = (!circle2d.IsNull());
2556
2557       if ( !isLine && !isCirc) // Check if the EDGE is close to a line
2558       {
2559         Bnd_B2d bndBox;
2560         SMDS_NodeIteratorPtr nIt = smDS->GetNodes();
2561         while ( nIt->more() )
2562           bndBox.Add( helper.GetNodeUV( F, nIt->next() ));
2563         gp_XY size = bndBox.CornerMax() - bndBox.CornerMin();
2564
2565         const double lineTol = 1e-2 * sqrt( bndBox.SquareExtent() );
2566         for ( int i = 0; i < 2 && !isLine; ++i )
2567           isLine = ( size.Coord( i+1 ) <= lineTol );
2568       }
2569       if ( !isLine && !isCirc && iTo-iFrom > 2) // Check if the EDGE is close to a circle
2570       {
2571         // TODO
2572       }
2573       if ( isLine )
2574       {
2575         line = new Geom_Line( gp::OX() ); // only type does matter
2576       }
2577       else if ( isCirc )
2578       {
2579         gp_Pnt2d p = circle2d->Location();
2580         gp_Ax2 ax( gp_Pnt( p.X(), p.Y(), 0), gp::DX());
2581         circle = new Geom_Circle( ax, 1.); // only center position does matter
2582       }
2583     }
2584
2585     Handle(Geom_Curve)& res = _edge2curve[ eIndex ];
2586     if ( isLine )
2587       res = line;
2588     else if ( isCirc )
2589       res = circle;
2590
2591     return res;
2592   }
2593   return i2curve->second;
2594 }
2595
2596 //================================================================================
2597 /*!
2598  * \brief smooth _LayerEdge's on a staight EDGE or circular EDGE
2599  */
2600 //================================================================================
2601
2602 bool _ViscousBuilder::smoothAnalyticEdge( _SolidData&           data,
2603                                           const int             iFrom,
2604                                           const int             iTo,
2605                                           Handle(Geom_Surface)& surface,
2606                                           const TopoDS_Face&    F,
2607                                           SMESH_MesherHelper&   helper)
2608 {
2609   TopoDS_Shape S = helper.GetSubShapeByNode( data._edges[ iFrom ]->_nodes[0],
2610                                              helper.GetMeshDS());
2611   TopoDS_Edge E = TopoDS::Edge( S );
2612
2613   Handle(Geom_Curve) curve = data.CurveForSmooth( E, iFrom, iTo, surface, F, helper );
2614   if ( curve.IsNull() ) return false;
2615
2616   // compute a relative length of segments
2617   vector< double > len( iTo-iFrom+1 );
2618   {
2619     double curLen, prevLen = len[0] = 1.0;
2620     for ( int i = iFrom; i < iTo; ++i )
2621     {
2622       curLen = prevLen * data._edges[i]->_2neibors->_wgt[0] / data._edges[i]->_2neibors->_wgt[1];
2623       len[i-iFrom+1] = len[i-iFrom] + curLen;
2624       prevLen = curLen;
2625     }
2626   }
2627
2628   if ( curve->IsKind( STANDARD_TYPE( Geom_Line )))
2629   {
2630     if ( F.IsNull() ) // 3D
2631     {
2632       SMESH_TNodeXYZ p0( data._edges[iFrom]->_2neibors->_nodes[0]);
2633       SMESH_TNodeXYZ p1( data._edges[iTo-1]->_2neibors->_nodes[1]);
2634       for ( int i = iFrom; i < iTo; ++i )
2635       {
2636         double r = len[i-iFrom] / len.back();
2637         gp_XYZ newPos = p0 * ( 1. - r ) + p1 * r;
2638         data._edges[i]->_pos.back() = newPos;
2639         SMDS_MeshNode* tgtNode = const_cast<SMDS_MeshNode*>( data._edges[i]->_nodes.back() );
2640         tgtNode->setXYZ( newPos.X(), newPos.Y(), newPos.Z() );
2641         dumpMove( tgtNode );
2642       }
2643     }
2644     else
2645     {
2646       gp_XY uv0 = helper.GetNodeUV( F, data._edges[iFrom]->_2neibors->_nodes[0]);
2647       gp_XY uv1 = helper.GetNodeUV( F, data._edges[iTo-1]->_2neibors->_nodes[1]);
2648       if ( data._edges[iFrom]->_2neibors->_nodes[0] ==
2649            data._edges[iTo-1]->_2neibors->_nodes[1] ) // closed edge
2650       {
2651         int iPeriodic = helper.GetPeriodicIndex();
2652         if ( iPeriodic == 1 || iPeriodic == 2 )
2653         {
2654           uv1.SetCoord( iPeriodic, helper.GetOtherParam( uv1.Coord( iPeriodic )));
2655           if ( uv0.Coord( iPeriodic ) > uv1.Coord( iPeriodic ))
2656             std::swap( uv0, uv1 );
2657         }
2658       }
2659       const gp_XY rangeUV = uv1 - uv0;
2660       for ( int i = iFrom; i < iTo; ++i )
2661       {
2662         double r = len[i-iFrom] / len.back();
2663         gp_XY newUV = uv0 + r * rangeUV;
2664         data._edges[i]->_pos.back().SetCoord( newUV.X(), newUV.Y(), 0 );
2665
2666         gp_Pnt newPos = surface->Value( newUV.X(), newUV.Y() );
2667         SMDS_MeshNode* tgtNode = const_cast<SMDS_MeshNode*>( data._edges[i]->_nodes.back() );
2668         tgtNode->setXYZ( newPos.X(), newPos.Y(), newPos.Z() );
2669         dumpMove( tgtNode );
2670
2671         SMDS_FacePosition* pos = static_cast<SMDS_FacePosition*>( tgtNode->GetPosition() );
2672         pos->SetUParameter( newUV.X() );
2673         pos->SetVParameter( newUV.Y() );
2674       }
2675     }
2676     return true;
2677   }
2678
2679   if ( curve->IsKind( STANDARD_TYPE( Geom_Circle )))
2680   {
2681     Handle(Geom_Circle) circle = Handle(Geom_Circle)::DownCast( curve );
2682     gp_Pnt center3D = circle->Location();
2683
2684     if ( F.IsNull() ) // 3D
2685     {
2686       if ( data._edges[iFrom]->_2neibors->_nodes[0] ==
2687            data._edges[iTo-1]->_2neibors->_nodes[1] )
2688         return true; // closed EDGE - nothing to do
2689
2690       return false; // TODO ???
2691     }
2692     else // 2D
2693     {
2694       const gp_XY center( center3D.X(), center3D.Y() );
2695
2696       gp_XY uv0 = helper.GetNodeUV( F, data._edges[iFrom]->_2neibors->_nodes[0]);
2697       gp_XY uvM = helper.GetNodeUV( F, data._edges[iFrom]->_nodes.back());
2698       gp_XY uv1 = helper.GetNodeUV( F, data._edges[iTo-1]->_2neibors->_nodes[1]);
2699       gp_Vec2d vec0( center, uv0 );
2700       gp_Vec2d vecM( center, uvM );
2701       gp_Vec2d vec1( center, uv1 );
2702       double uLast = vec0.Angle( vec1 ); // -PI - +PI
2703       double uMidl = vec0.Angle( vecM );
2704       if ( uLast * uMidl < 0. )
2705         uLast += ( uMidl > 0 ? +2. : -2. ) * M_PI;
2706       const double radius = 0.5 * ( vec0.Magnitude() + vec1.Magnitude() );
2707
2708       gp_Ax2d   axis( center, vec0 );
2709       gp_Circ2d circ( axis, radius );
2710       for ( int i = iFrom; i < iTo; ++i )
2711       {
2712         double    newU = uLast * len[i-iFrom] / len.back();
2713         gp_Pnt2d newUV = ElCLib::Value( newU, circ );
2714         data._edges[i]->_pos.back().SetCoord( newUV.X(), newUV.Y(), 0 );
2715
2716         gp_Pnt newPos = surface->Value( newUV.X(), newUV.Y() );
2717         SMDS_MeshNode* tgtNode = const_cast<SMDS_MeshNode*>( data._edges[i]->_nodes.back() );
2718         tgtNode->setXYZ( newPos.X(), newPos.Y(), newPos.Z() );
2719         dumpMove( tgtNode );
2720
2721         SMDS_FacePosition* pos = static_cast<SMDS_FacePosition*>( tgtNode->GetPosition() );
2722         pos->SetUParameter( newUV.X() );
2723         pos->SetVParameter( newUV.Y() );
2724       }
2725     }
2726     return true;
2727   }
2728
2729   return false;
2730 }
2731
2732 //================================================================================
2733 /*!
2734  * \brief Modify normals of _LayerEdge's on EDGE's to avoid intersection with
2735  * _LayerEdge's on neighbor EDGE's
2736  */
2737 //================================================================================
2738
2739 bool _ViscousBuilder::updateNormals( _SolidData&         data,
2740                                      SMESH_MesherHelper& helper )
2741 {
2742   // make temporary quadrangles got by extrusion of
2743   // mesh edges along _LayerEdge._normal's
2744
2745   vector< const SMDS_MeshElement* > tmpFaces;
2746   {
2747     set< SMESH_TLink > extrudedLinks; // contains target nodes
2748     vector< const SMDS_MeshNode*> nodes(4); // of a tmp mesh face
2749
2750     dumpFunction(SMESH_Comment("makeTmpFacesOnEdges")<<data._index);
2751     for ( size_t i = 0; i < data._edges.size(); ++i )
2752     {
2753       _LayerEdge* edge = data._edges[i];
2754       if ( !edge->IsOnEdge() || !edge->_sWOL.IsNull() ) continue;
2755       const SMDS_MeshNode* tgt1 = edge->_nodes.back();
2756       for ( int j = 0; j < 2; ++j ) // loop on _2NearEdges
2757       {
2758         const SMDS_MeshNode* tgt2 = edge->_2neibors->_nodes[j];
2759         pair< set< SMESH_TLink >::iterator, bool > link_isnew =
2760           extrudedLinks.insert( SMESH_TLink( tgt1, tgt2 ));
2761         if ( !link_isnew.second )
2762         {
2763           extrudedLinks.erase( link_isnew.first );
2764           continue; // already extruded and will no more encounter
2765         }
2766         // look for a _LayerEdge containg tgt2
2767 //         _LayerEdge* neiborEdge = 0;
2768 //         size_t di = 0; // check _edges[i+di] and _edges[i-di]
2769 //         while ( !neiborEdge && ++di <= data._edges.size() )
2770 //         {
2771 //           if ( i+di < data._edges.size() && data._edges[i+di]->_nodes.back() == tgt2 )
2772 //             neiborEdge = data._edges[i+di];
2773 //           else if ( di <= i && data._edges[i-di]->_nodes.back() == tgt2 )
2774 //             neiborEdge = data._edges[i-di];
2775 //         }
2776 //         if ( !neiborEdge )
2777 //           return error("updateNormals(): neighbor _LayerEdge not found", data._index);
2778         _LayerEdge* neiborEdge = edge->_2neibors->_edges[j];
2779
2780         TmpMeshFaceOnEdge* f = new TmpMeshFaceOnEdge( edge, neiborEdge, --_tmpFaceID );
2781         tmpFaces.push_back( f );
2782
2783         dumpCmd(SMESH_Comment("mesh.AddFace([ ")
2784                 <<f->_nn[0]->GetID()<<", "<<f->_nn[1]->GetID()<<", "
2785                 <<f->_nn[2]->GetID()<<", "<<f->_nn[3]->GetID()<<" ])");
2786       }
2787     }
2788     dumpFunctionEnd();
2789   }
2790   // Check if _LayerEdge's based on EDGE's intersects tmpFaces.
2791   // Perform two loops on _LayerEdge on EDGE's:
2792   // 1) to find and fix intersection
2793   // 2) to check that no new intersection appears as result of 1)
2794
2795   SMDS_ElemIteratorPtr fIt( new SMDS_ElementVectorIterator( tmpFaces.begin(),
2796                                                             tmpFaces.end()));
2797   auto_ptr<SMESH_ElementSearcher> searcher
2798     ( SMESH_MeshAlgos::GetElementSearcher( *getMeshDS(), fIt ));
2799
2800   // 1) Find intersections
2801   double dist;
2802   const SMDS_MeshElement* face;
2803   typedef map< _LayerEdge*, set< _LayerEdge*, _LayerEdgeCmp >, _LayerEdgeCmp > TLEdge2LEdgeSet;
2804   TLEdge2LEdgeSet edge2CloseEdge;
2805
2806   const double eps = data._epsilon * data._epsilon;
2807   for ( size_t i = 0; i < data._edges.size(); ++i )
2808   {
2809     _LayerEdge* edge = data._edges[i];
2810     if ( !edge->IsOnEdge() || !edge->_sWOL.IsNull() ) continue;
2811     if ( edge->FindIntersection( *searcher, dist, eps, &face ))
2812     {
2813       const TmpMeshFaceOnEdge* f = (const TmpMeshFaceOnEdge*) face;
2814       set< _LayerEdge*, _LayerEdgeCmp > & ee = edge2CloseEdge[ edge ];
2815       ee.insert( f->_le1 );
2816       ee.insert( f->_le2 );
2817       if ( f->_le1->IsOnEdge() && f->_le1->_sWOL.IsNull() ) 
2818         edge2CloseEdge[ f->_le1 ].insert( edge );
2819       if ( f->_le2->IsOnEdge() && f->_le2->_sWOL.IsNull() ) 
2820         edge2CloseEdge[ f->_le2 ].insert( edge );
2821     }
2822   }
2823
2824   // Set _LayerEdge._normal
2825
2826   if ( !edge2CloseEdge.empty() )
2827   {
2828     dumpFunction(SMESH_Comment("updateNormals")<<data._index);
2829
2830     TLEdge2LEdgeSet::iterator e2ee = edge2CloseEdge.begin();
2831     for ( ; e2ee != edge2CloseEdge.end(); ++e2ee )
2832     {
2833       _LayerEdge* edge1       = e2ee->first;
2834       _LayerEdge* edge2       = 0;
2835       set< _LayerEdge*, _LayerEdgeCmp >& ee  = e2ee->second;
2836
2837       // find EDGEs the edges reside
2838       TopoDS_Edge E1, E2;
2839       TopoDS_Shape S = helper.GetSubShapeByNode( edge1->_nodes[0], getMeshDS() );
2840       if ( S.ShapeType() != TopAbs_EDGE )
2841         continue; // TODO: find EDGE by VERTEX
2842       E1 = TopoDS::Edge( S );
2843       set< _LayerEdge*, _LayerEdgeCmp >::iterator eIt = ee.begin();
2844       while ( E2.IsNull() && eIt != ee.end())
2845       {
2846         _LayerEdge* e2 = *eIt++;
2847         TopoDS_Shape S = helper.GetSubShapeByNode( e2->_nodes[0], getMeshDS() );
2848         if ( S.ShapeType() == TopAbs_EDGE )
2849           E2 = TopoDS::Edge( S ), edge2 = e2;
2850       }
2851       if ( E2.IsNull() ) continue; // TODO: find EDGE by VERTEX
2852
2853       // find 3 FACEs sharing 2 EDGEs
2854
2855       TopoDS_Face FF1[2], FF2[2];
2856       PShapeIteratorPtr fIt = helper.GetAncestors(E1, *_mesh, TopAbs_FACE);
2857       while ( fIt->more() && FF1[1].IsNull())
2858       {
2859         const TopoDS_Face *F = (const TopoDS_Face*) fIt->next();
2860         if ( helper.IsSubShape( *F, data._solid))
2861           FF1[ FF1[0].IsNull() ? 0 : 1 ] = *F;
2862       }
2863       fIt = helper.GetAncestors(E2, *_mesh, TopAbs_FACE);
2864       while ( fIt->more() && FF2[1].IsNull())
2865       {
2866         const TopoDS_Face *F = (const TopoDS_Face*) fIt->next();
2867         if ( helper.IsSubShape( *F, data._solid))
2868           FF2[ FF2[0].IsNull() ? 0 : 1 ] = *F;
2869       }
2870       // exclude a FACE common to E1 and E2 (put it at [1] in FF* )
2871       if ( FF1[0].IsSame( FF2[0]) || FF1[0].IsSame( FF2[1]))
2872         std::swap( FF1[0], FF1[1] );
2873       if ( FF2[0].IsSame( FF1[0]) )
2874         std::swap( FF2[0], FF2[1] );
2875       if ( FF1[0].IsNull() || FF2[0].IsNull() )
2876         continue;
2877
2878 //       // get a new normal for edge1
2879       bool ok;
2880       gp_Vec dir1 = edge1->_normal, dir2 = edge2->_normal;
2881       if ( edge1->_cosin < 0 )
2882         dir1 = getFaceDir( FF1[0], E1, edge1->_nodes[0], helper, ok ).Normalized();
2883       if ( edge2->_cosin < 0 )
2884         dir2 = getFaceDir( FF2[0], E2, edge2->_nodes[0], helper, ok ).Normalized();
2885       //      gp_Vec dir1 = getFaceDir( FF1[0], E1, edge1->_nodes[0], helper, ok );
2886 //       gp_Vec dir2 = getFaceDir( FF2[0], E2, edge2->_nodes[0], helper, ok2 );
2887 //       double wgt1 = ( edge1->_cosin + 1 ) / ( edge1->_cosin + edge2->_cosin + 2 );
2888 //       double wgt2 = ( edge2->_cosin + 1 ) / ( edge1->_cosin + edge2->_cosin + 2 );
2889 //       gp_Vec newNorm = wgt1 * dir1 + wgt2 * dir2;
2890 //       newNorm.Normalize();
2891
2892       double wgt1 = ( edge1->_cosin + 1 ) / ( edge1->_cosin + edge2->_cosin + 2 );
2893       double wgt2 = ( edge2->_cosin + 1 ) / ( edge1->_cosin + edge2->_cosin + 2 );
2894       gp_Vec newNorm = wgt1 * dir1 + wgt2 * dir2;
2895       newNorm.Normalize();
2896
2897       edge1->_normal = newNorm.XYZ();
2898
2899       // update data of edge1 depending on _normal
2900       const SMDS_MeshNode *n1, *n2;
2901       n1 = edge1->_2neibors->_edges[0]->_nodes[0];
2902       n2 = edge1->_2neibors->_edges[1]->_nodes[0];
2903       //if ( !findNeiborsOnEdge( edge1, n1, n2, data ))
2904       //continue;
2905       edge1->SetDataByNeighbors( n1, n2, helper );
2906       gp_Vec dirInFace;
2907       if ( edge1->_cosin < 0 )
2908         dirInFace = dir1;
2909       else
2910         getFaceDir( FF1[0], E1, edge1->_nodes[0], helper, ok );
2911       double angle = dir1.Angle( edge1->_normal ); // [0,PI]
2912       edge1->SetCosin( cos( angle ));
2913
2914       // limit data._stepSize
2915       if ( edge1->_cosin > 0.1 )
2916       {
2917         SMDS_ElemIteratorPtr fIt = edge1->_nodes[0]->GetInverseElementIterator(SMDSAbs_Face);
2918         while ( fIt->more() )
2919           limitStepSize( data, fIt->next(), edge1->_cosin );
2920       }
2921       // set new XYZ of target node
2922       edge1->InvalidateStep( 1 );
2923       edge1->_len = 0;
2924       edge1->SetNewLength( data._stepSize, helper );
2925     }
2926
2927     // Update normals and other dependent data of not intersecting _LayerEdge's
2928     // neighboring the intersecting ones
2929
2930     for ( e2ee = edge2CloseEdge.begin(); e2ee != edge2CloseEdge.end(); ++e2ee )
2931     {
2932       _LayerEdge* edge1 = e2ee->first;
2933       if ( !edge1->_2neibors )
2934         continue;
2935       for ( int j = 0; j < 2; ++j ) // loop on 2 neighbors
2936       {
2937         _LayerEdge* neighbor = edge1->_2neibors->_edges[j];
2938         if ( edge2CloseEdge.count ( neighbor ))
2939           continue; // j-th neighbor is also intersected
2940         _LayerEdge* prevEdge = edge1;
2941         const int nbSteps = 6;
2942         for ( int step = nbSteps; step; --step ) // step from edge1 in j-th direction
2943         {
2944           if ( !neighbor->_2neibors )
2945             break; // neighbor is on VERTEX
2946           int iNext = 0;
2947           _LayerEdge* nextEdge = neighbor->_2neibors->_edges[iNext];
2948           if ( nextEdge == prevEdge )
2949             nextEdge = neighbor->_2neibors->_edges[ ++iNext ];
2950 //           const double&  wgtPrev = neighbor->_2neibors->_wgt[1-iNext];
2951 //           const double&  wgtNext = neighbor->_2neibors->_wgt[iNext];
2952           double r = double(step-1)/nbSteps;
2953           if ( !nextEdge->_2neibors )
2954             r = 0.5;
2955
2956           gp_XYZ newNorm = prevEdge->_normal * r + nextEdge->_normal * (1-r);
2957           newNorm.Normalize();
2958
2959           neighbor->_normal = newNorm;
2960           neighbor->SetCosin( prevEdge->_cosin * r + nextEdge->_cosin * (1-r) );
2961           neighbor->SetDataByNeighbors( prevEdge->_nodes[0], nextEdge->_nodes[0], helper );
2962
2963           neighbor->InvalidateStep( 1 );
2964           neighbor->_len = 0;
2965           neighbor->SetNewLength( data._stepSize, helper );
2966
2967           // goto the next neighbor
2968           prevEdge = neighbor;
2969           neighbor = nextEdge;
2970         }
2971       }
2972     }
2973     dumpFunctionEnd();
2974   }
2975   // 2) Check absence of intersections
2976   // TODO?
2977
2978   for ( size_t i = 0 ; i < tmpFaces.size(); ++i )
2979     delete tmpFaces[i];
2980
2981   return true;
2982 }
2983
2984 //================================================================================
2985 /*!
2986  * \brief Looks for intersection of it's last segment with faces
2987  *  \param distance - returns shortest distance from the last node to intersection
2988  */
2989 //================================================================================
2990
2991 bool _LayerEdge::FindIntersection( SMESH_ElementSearcher&   searcher,
2992                                    double &                 distance,
2993                                    const double&            epsilon,
2994                                    const SMDS_MeshElement** face)
2995 {
2996   vector< const SMDS_MeshElement* > suspectFaces;
2997   double segLen;
2998   gp_Ax1 lastSegment = LastSegment(segLen);
2999   searcher.GetElementsNearLine( lastSegment, SMDSAbs_Face, suspectFaces );
3000
3001   bool segmentIntersected = false;
3002   distance = Precision::Infinite();
3003   int iFace = -1; // intersected face
3004   for ( size_t j = 0 ; j < suspectFaces.size() && !segmentIntersected; ++j )
3005   {
3006     const SMDS_MeshElement* face = suspectFaces[j];
3007     if ( face->GetNodeIndex( _nodes.back() ) >= 0 ||
3008          face->GetNodeIndex( _nodes[0]     ) >= 0 )
3009       continue; // face sharing _LayerEdge node
3010     const int nbNodes = face->NbCornerNodes();
3011     bool intFound = false;
3012     double dist;
3013     SMDS_MeshElement::iterator nIt = face->begin_nodes();
3014     if ( nbNodes == 3 )
3015     {
3016       intFound = SegTriaInter( lastSegment, *nIt++, *nIt++, *nIt++, dist, epsilon );
3017     }
3018     else
3019     {
3020       const SMDS_MeshNode* tria[3];
3021       tria[0] = *nIt++;
3022       tria[1] = *nIt++;;
3023       for ( int n2 = 2; n2 < nbNodes && !intFound; ++n2 )
3024       {
3025         tria[2] = *nIt++;
3026         intFound = SegTriaInter(lastSegment, tria[0], tria[1], tria[2], dist, epsilon );
3027         tria[1] = tria[2];
3028       }
3029     }
3030     if ( intFound )
3031     {
3032       if ( dist < segLen*(1.01) && dist > -(_len-segLen) )
3033         segmentIntersected = true;
3034       if ( distance > dist )
3035         distance = dist, iFace = j;
3036     }
3037   }
3038   if ( iFace != -1 && face ) *face = suspectFaces[iFace];
3039 //   if ( distance && iFace > -1 )
3040 //   {
3041 //     // distance is used to limit size of inflation step which depends on
3042 //     // whether the intersected face bears viscous layers or not
3043 //     bool faceHasVL = suspectFaces[iFace]->GetID() < 1;
3044 //     if ( faceHasVL )
3045 //       *distance /= 2;
3046 //   }
3047   if ( segmentIntersected )
3048   {
3049 #ifdef __myDEBUG
3050     SMDS_MeshElement::iterator nIt = suspectFaces[iFace]->begin_nodes();
3051     gp_XYZ intP( lastSegment.Location().XYZ() + lastSegment.Direction().XYZ() * distance );
3052     cout << "nodes: tgt " << _nodes.back()->GetID() << " src " << _nodes[0]->GetID()
3053          << ", intersection with face ("
3054          << (*nIt++)->GetID()<<" "<< (*nIt++)->GetID()<<" "<< (*nIt++)->GetID()
3055          << ") at point (" << intP.X() << ", " << intP.Y() << ", " << intP.Z()
3056          << ") distance = " << distance - segLen<< endl;
3057 #endif
3058   }
3059
3060   distance -= segLen;
3061
3062   return segmentIntersected;
3063 }
3064
3065 //================================================================================
3066 /*!
3067  * \brief Returns size and direction of the last segment
3068  */
3069 //================================================================================
3070
3071 gp_Ax1 _LayerEdge::LastSegment(double& segLen) const
3072 {
3073   // find two non-coincident positions
3074   gp_XYZ orig = _pos.back();
3075   gp_XYZ dir;
3076   int iPrev = _pos.size() - 2;
3077   while ( iPrev >= 0 )
3078   {
3079     dir = orig - _pos[iPrev];
3080     if ( dir.SquareModulus() > 1e-100 )
3081       break;
3082     else
3083       iPrev--;
3084   }
3085
3086   // make gp_Ax1
3087   gp_Ax1 segDir;
3088   if ( iPrev < 0 )
3089   {
3090     segDir.SetLocation( SMESH_TNodeXYZ( _nodes[0] ));
3091     segDir.SetDirection( _normal );
3092     segLen = 0;
3093   }
3094   else
3095   {
3096     gp_Pnt pPrev = _pos[ iPrev ];
3097     if ( !_sWOL.IsNull() )
3098     {
3099       TopLoc_Location loc;
3100       if ( _sWOL.ShapeType() == TopAbs_EDGE )
3101       {
3102         double f,l;
3103         Handle(Geom_Curve) curve = BRep_Tool::Curve( TopoDS::Edge( _sWOL ), loc, f,l);
3104         pPrev = curve->Value( pPrev.X() ).Transformed( loc );
3105       }
3106       else
3107       {
3108         Handle(Geom_Surface) surface = BRep_Tool::Surface( TopoDS::Face(_sWOL), loc );
3109         pPrev = surface->Value( pPrev.X(), pPrev.Y() ).Transformed( loc );
3110       }
3111       dir = SMESH_TNodeXYZ( _nodes.back() ) - pPrev.XYZ();
3112     }
3113     segDir.SetLocation( pPrev );
3114     segDir.SetDirection( dir );
3115     segLen = dir.Modulus();
3116   }
3117
3118   return segDir;
3119 }
3120
3121 //================================================================================
3122 /*!
3123  * \brief Test intersection of the last segment with a given triangle
3124  *   using Moller-Trumbore algorithm
3125  * Intersection is detected if distance to intersection is less than _LayerEdge._len
3126  */
3127 //================================================================================
3128
3129 bool _LayerEdge::SegTriaInter( const gp_Ax1&        lastSegment,
3130                                const SMDS_MeshNode* n0,
3131                                const SMDS_MeshNode* n1,
3132                                const SMDS_MeshNode* n2,
3133                                double&              t,
3134                                const double&        EPSILON) const
3135 {
3136   //const double EPSILON = 1e-6;
3137
3138   gp_XYZ orig = lastSegment.Location().XYZ();
3139   gp_XYZ dir  = lastSegment.Direction().XYZ();
3140
3141   SMESH_TNodeXYZ vert0( n0 );
3142   SMESH_TNodeXYZ vert1( n1 );
3143   SMESH_TNodeXYZ vert2( n2 );
3144
3145   /* calculate distance from vert0 to ray origin */
3146   gp_XYZ tvec = orig - vert0;
3147
3148   if ( tvec * dir > EPSILON )
3149     // intersected face is at back side of the temporary face this _LayerEdge belongs to
3150     return false;
3151
3152   gp_XYZ edge1 = vert1 - vert0;
3153   gp_XYZ edge2 = vert2 - vert0;
3154
3155   /* begin calculating determinant - also used to calculate U parameter */
3156   gp_XYZ pvec = dir ^ edge2;
3157
3158   /* if determinant is near zero, ray lies in plane of triangle */
3159   double det = edge1 * pvec;
3160
3161   if (det > -EPSILON && det < EPSILON)
3162     return 0;
3163   double inv_det = 1.0 / det;
3164
3165   /* calculate U parameter and test bounds */
3166   double u = ( tvec * pvec ) * inv_det;
3167   if (u < 0.0 || u > 1.0)
3168     return 0;
3169
3170   /* prepare to test V parameter */
3171   gp_XYZ qvec = tvec ^ edge1;
3172
3173   /* calculate V parameter and test bounds */
3174   double v = (dir * qvec) * inv_det;
3175   if ( v < 0.0 || u + v > 1.0 )
3176     return 0;
3177
3178   /* calculate t, ray intersects triangle */
3179   t = (edge2 * qvec) * inv_det;
3180
3181   //   if (det < EPSILON)
3182   //     return false;
3183
3184   //   /* calculate distance from vert0 to ray origin */
3185   //   gp_XYZ tvec = orig - vert0;
3186
3187   //   /* calculate U parameter and test bounds */
3188   //   double u = tvec * pvec;
3189   //   if (u < 0.0 || u > det)
3190 //     return 0;
3191
3192 //   /* prepare to test V parameter */
3193 //   gp_XYZ qvec = tvec ^ edge1;
3194
3195 //   /* calculate V parameter and test bounds */
3196 //   double v = dir * qvec;
3197 //   if (v < 0.0 || u + v > det)
3198 //     return 0;
3199
3200 //   /* calculate t, scale parameters, ray intersects triangle */
3201 //   double t = edge2 * qvec;
3202 //   double inv_det = 1.0 / det;
3203 //   t *= inv_det;
3204 //   //u *= inv_det;
3205 //   //v *= inv_det;
3206
3207   return true;
3208 }
3209
3210 //================================================================================
3211 /*!
3212  * \brief Perform smooth of _LayerEdge's based on EDGE's
3213  *  \retval bool - true if node has been moved
3214  */
3215 //================================================================================
3216
3217 bool _LayerEdge::SmoothOnEdge(Handle(Geom_Surface)& surface,
3218                               const TopoDS_Face&    F,
3219                               SMESH_MesherHelper&   helper)
3220 {
3221   ASSERT( IsOnEdge() );
3222
3223   SMDS_MeshNode* tgtNode = const_cast<SMDS_MeshNode*>( _nodes.back() );
3224   SMESH_TNodeXYZ oldPos( tgtNode );
3225   double dist01, distNewOld;
3226   
3227   SMESH_TNodeXYZ p0( _2neibors->_nodes[0]);
3228   SMESH_TNodeXYZ p1( _2neibors->_nodes[1]);
3229   dist01 = p0.Distance( _2neibors->_nodes[1] );
3230
3231   gp_Pnt newPos = p0 * _2neibors->_wgt[0] + p1 * _2neibors->_wgt[1];
3232   double lenDelta = 0;
3233   if ( _curvature )
3234   {
3235     //lenDelta = _curvature->lenDelta( _len );
3236     lenDelta = _curvature->lenDeltaByDist( dist01 );
3237     newPos.ChangeCoord() += _normal * lenDelta;
3238   }
3239
3240   distNewOld = newPos.Distance( oldPos );
3241
3242   if ( F.IsNull() )
3243   {
3244     if ( _2neibors->_plnNorm )
3245     {
3246       // put newPos on the plane defined by source node and _plnNorm
3247       gp_XYZ new2src = SMESH_TNodeXYZ( _nodes[0] ) - newPos.XYZ();
3248       double new2srcProj = (*_2neibors->_plnNorm) * new2src;
3249       newPos.ChangeCoord() += (*_2neibors->_plnNorm) * new2srcProj;
3250     }
3251     tgtNode->setXYZ( newPos.X(), newPos.Y(), newPos.Z() );
3252     _pos.back() = newPos.XYZ();
3253   }
3254   else
3255   {
3256     tgtNode->setXYZ( newPos.X(), newPos.Y(), newPos.Z() );
3257     gp_XY uv( Precision::Infinite(), 0 );
3258     helper.CheckNodeUV( F, tgtNode, uv, 1e-10, /*force=*/true );
3259     _pos.back().SetCoord( uv.X(), uv.Y(), 0 );
3260
3261     newPos = surface->Value( uv.X(), uv.Y() );
3262     tgtNode->setXYZ( newPos.X(), newPos.Y(), newPos.Z() );
3263   }
3264
3265   if ( _curvature && lenDelta < 0 )
3266   {
3267     gp_Pnt prevPos( _pos[ _pos.size()-2 ]);
3268     _len -= prevPos.Distance( oldPos );
3269     _len += prevPos.Distance( newPos );
3270   }
3271   bool moved = distNewOld > dist01/50;
3272   //if ( moved )
3273   dumpMove( tgtNode ); // debug
3274
3275   return moved;
3276 }
3277
3278 //================================================================================
3279 /*!
3280  * \brief Perform laplacian smooth in 3D of nodes inflated from FACE
3281  *  \retval bool - true if _tgtNode has been moved
3282  */
3283 //================================================================================
3284
3285 bool _LayerEdge::Smooth(int& badNb)
3286 {
3287   if ( _simplices.size() < 2 )
3288     return false; // _LayerEdge inflated along EDGE or FACE
3289
3290   // compute new position for the last _pos
3291   gp_XYZ newPos (0,0,0);
3292   for ( size_t i = 0; i < _simplices.size(); ++i )
3293     newPos += SMESH_TNodeXYZ( _simplices[i]._nPrev );
3294   newPos /= _simplices.size();
3295
3296   if ( _curvature )
3297     newPos += _normal * _curvature->lenDelta( _len );
3298
3299   gp_Pnt prevPos( _pos[ _pos.size()-2 ]);
3300 //   if ( _cosin < -0.1)
3301 //   {
3302 //     // Avoid decreasing length of edge on concave surface
3303 //     //gp_Vec oldMove( _pos[ _pos.size()-2 ], _pos.back() );
3304 //     gp_Vec newMove( prevPos, newPos );
3305 //     newPos = _pos.back() + newMove.XYZ();
3306 //   }
3307 //   else if ( _cosin > 0.3 )
3308 //   {
3309 //     // Avoid increasing length of edge too much
3310
3311 //   }
3312   // count quality metrics (orientation) of tetras around _tgtNode
3313   int nbOkBefore = 0;
3314   SMESH_TNodeXYZ tgtXYZ( _nodes.back() );
3315   for ( size_t i = 0; i < _simplices.size(); ++i )
3316     nbOkBefore += _simplices[i].IsForward( _nodes[0], &tgtXYZ );
3317
3318   int nbOkAfter = 0;
3319   for ( size_t i = 0; i < _simplices.size(); ++i )
3320     nbOkAfter += _simplices[i].IsForward( _nodes[0], &newPos );
3321
3322   if ( nbOkAfter < nbOkBefore )
3323     return false;
3324
3325   SMDS_MeshNode* n = const_cast< SMDS_MeshNode* >( _nodes.back() );
3326
3327   _len -= prevPos.Distance(SMESH_TNodeXYZ( n ));
3328   _len += prevPos.Distance(newPos);
3329
3330   n->setXYZ( newPos.X(), newPos.Y(), newPos.Z());
3331   _pos.back() = newPos;
3332
3333   badNb += _simplices.size() - nbOkAfter;
3334
3335   dumpMove( n );
3336
3337   return true;
3338 }
3339
3340 //================================================================================
3341 /*!
3342  * \brief Add a new segment to _LayerEdge during inflation
3343  */
3344 //================================================================================
3345
3346 void _LayerEdge::SetNewLength( double len, SMESH_MesherHelper& helper )
3347 {
3348   if ( _len - len > -1e-6 )
3349   {
3350     _pos.push_back( _pos.back() );
3351     return;
3352   }
3353
3354   SMDS_MeshNode* n = const_cast< SMDS_MeshNode*>( _nodes.back() );
3355   SMESH_TNodeXYZ oldXYZ( n );
3356   gp_XYZ nXYZ = oldXYZ + _normal * ( len - _len ) * _lenFactor;
3357   n->setXYZ( nXYZ.X(), nXYZ.Y(), nXYZ.Z() );
3358
3359   _pos.push_back( nXYZ );
3360   _len = len;
3361   if ( !_sWOL.IsNull() )
3362   {
3363     double distXYZ[4];
3364     if ( _sWOL.ShapeType() == TopAbs_EDGE )
3365     {
3366       double u = Precision::Infinite(); // to force projection w/o distance check
3367       helper.CheckNodeU( TopoDS::Edge( _sWOL ), n, u, 1e-10, /*force=*/true, distXYZ );
3368       _pos.back().SetCoord( u, 0, 0 );
3369       SMDS_EdgePosition* pos = static_cast<SMDS_EdgePosition*>( n->GetPosition() );
3370       pos->SetUParameter( u );
3371     }
3372     else //  TopAbs_FACE
3373     {
3374       gp_XY uv( Precision::Infinite(), 0 );
3375       helper.CheckNodeUV( TopoDS::Face( _sWOL ), n, uv, 1e-10, /*force=*/true, distXYZ );
3376       _pos.back().SetCoord( uv.X(), uv.Y(), 0 );
3377       SMDS_FacePosition* pos = static_cast<SMDS_FacePosition*>( n->GetPosition() );
3378       pos->SetUParameter( uv.X() );
3379       pos->SetVParameter( uv.Y() );
3380     }
3381     n->setXYZ( distXYZ[1], distXYZ[2], distXYZ[3]);
3382   }
3383   dumpMove( n ); //debug
3384 }
3385
3386 //================================================================================
3387 /*!
3388  * \brief Remove last inflation step
3389  */
3390 //================================================================================
3391
3392 void _LayerEdge::InvalidateStep( int curStep )
3393 {
3394   if ( _pos.size() > curStep )
3395   {
3396     _pos.resize( curStep );
3397     gp_Pnt nXYZ = _pos.back();
3398     SMDS_MeshNode* n = const_cast< SMDS_MeshNode*>( _nodes.back() );
3399     if ( !_sWOL.IsNull() )
3400     {
3401       TopLoc_Location loc;
3402       if ( _sWOL.ShapeType() == TopAbs_EDGE )
3403       {
3404         SMDS_EdgePosition* pos = static_cast<SMDS_EdgePosition*>( n->GetPosition() );
3405         pos->SetUParameter( nXYZ.X() );
3406         double f,l;
3407         Handle(Geom_Curve) curve = BRep_Tool::Curve( TopoDS::Edge( _sWOL ), loc, f,l);
3408         nXYZ = curve->Value( nXYZ.X() ).Transformed( loc );
3409       }
3410       else
3411       {
3412         SMDS_FacePosition* pos = static_cast<SMDS_FacePosition*>( n->GetPosition() );
3413         pos->SetUParameter( nXYZ.X() );
3414         pos->SetVParameter( nXYZ.Y() );
3415         Handle(Geom_Surface) surface = BRep_Tool::Surface( TopoDS::Face(_sWOL), loc );
3416         nXYZ = surface->Value( nXYZ.X(), nXYZ.Y() ).Transformed( loc );
3417       }
3418     }
3419     n->setXYZ( nXYZ.X(), nXYZ.Y(), nXYZ.Z() );
3420     dumpMove( n );
3421   }
3422 }
3423
3424 //================================================================================
3425 /*!
3426  * \brief Create layers of prisms
3427  */
3428 //================================================================================
3429
3430 bool _ViscousBuilder::refine(_SolidData& data)
3431 {
3432   SMESH_MesherHelper helper( *_mesh );
3433   helper.SetSubShape( data._solid );
3434   helper.SetElementsOnShape(false);
3435
3436   Handle(Geom_Curve) curve;
3437   Handle(Geom_Surface) surface;
3438   TopoDS_Edge geomEdge;
3439   TopoDS_Face geomFace;
3440   TopLoc_Location loc;
3441   double f,l, u/*, distXYZ[4]*/;
3442   gp_XY uv;
3443   bool isOnEdge;
3444
3445   for ( size_t i = 0; i < data._edges.size(); ++i )
3446   {
3447     _LayerEdge& edge = *data._edges[i];
3448
3449     // get accumulated length of segments
3450     vector< double > segLen( edge._pos.size() );
3451     segLen[0] = 0.0;
3452     for ( size_t j = 1; j < edge._pos.size(); ++j )
3453       segLen[j] = segLen[j-1] + (edge._pos[j-1] - edge._pos[j] ).Modulus();
3454
3455     // allocate memory for new nodes if it is not yet refined
3456     const SMDS_MeshNode* tgtNode = edge._nodes.back();
3457     if ( edge._nodes.size() == 2 )
3458     {
3459       edge._nodes.resize( data._hyp->GetNumberLayers() + 1, 0 );
3460       edge._nodes[1] = 0;
3461       edge._nodes.back() = tgtNode;
3462     }
3463     if ( !edge._sWOL.IsNull() )
3464     {
3465       isOnEdge = ( edge._sWOL.ShapeType() == TopAbs_EDGE );
3466       // restore position of the last node
3467 //       gp_Pnt p;
3468       if ( isOnEdge )
3469       {
3470         geomEdge = TopoDS::Edge( edge._sWOL );
3471         curve = BRep_Tool::Curve( geomEdge, loc, f,l);
3472 //         double u = helper.GetNodeU( tgtNode );
3473 //         p = curve->Value( u );
3474       }
3475       else
3476       {
3477         geomFace = TopoDS::Face( edge._sWOL );
3478         surface = BRep_Tool::Surface( geomFace, loc );
3479 //         gp_XY uv = helper.GetNodeUV( tgtNode );
3480 //         p = surface->Value( uv.X(), uv.Y() );
3481       }
3482 //       p.Transform( loc );
3483 //       const_cast< SMDS_MeshNode* >( tgtNode )->setXYZ( p.X(), p.Y(), p.Z() );
3484     }
3485     // calculate height of the first layer
3486     double h0;
3487     const double T = segLen.back(); //data._hyp.GetTotalThickness();
3488     const double f = data._hyp->GetStretchFactor();
3489     const int    N = data._hyp->GetNumberLayers();
3490     const double fPowN = pow( f, N );
3491     if ( fPowN - 1 <= numeric_limits<double>::min() )
3492       h0 = T / N;
3493     else
3494       h0 = T * ( f - 1 )/( fPowN - 1 );
3495
3496     const double zeroLen = std::numeric_limits<double>::min();
3497
3498     // create intermediate nodes
3499     double hSum = 0, hi = h0/f;
3500     size_t iSeg = 1;
3501     for ( size_t iStep = 1; iStep < edge._nodes.size(); ++iStep )
3502     {
3503       // compute an intermediate position
3504       hi *= f;
3505       hSum += hi;
3506       while ( hSum > segLen[iSeg] && iSeg < segLen.size()-1)
3507         ++iSeg;
3508       int iPrevSeg = iSeg-1;
3509       while ( fabs( segLen[iPrevSeg] - segLen[iSeg]) <= zeroLen && iPrevSeg > 0 )
3510         --iPrevSeg;
3511       double r = ( segLen[iSeg] - hSum ) / ( segLen[iSeg] - segLen[iPrevSeg] );
3512       gp_Pnt pos = r * edge._pos[iPrevSeg] + (1-r) * edge._pos[iSeg];
3513
3514       SMDS_MeshNode*& node = const_cast< SMDS_MeshNode*& >(edge._nodes[ iStep ]);
3515       if ( !edge._sWOL.IsNull() )
3516       {
3517         // compute XYZ by parameters <pos>
3518         if ( isOnEdge )
3519         {
3520           u = pos.X();
3521           pos = curve->Value( u ).Transformed(loc);
3522         }
3523         else
3524         {
3525           uv.SetCoord( pos.X(), pos.Y() );
3526           pos = surface->Value( pos.X(), pos.Y() ).Transformed(loc);
3527         }
3528       }
3529       // create or update the node
3530       if ( !node )
3531       {
3532         node = helper.AddNode( pos.X(), pos.Y(), pos.Z());
3533         if ( !edge._sWOL.IsNull() )
3534         {
3535           if ( isOnEdge )
3536             getMeshDS()->SetNodeOnEdge( node, geomEdge, u );
3537           else
3538             getMeshDS()->SetNodeOnFace( node, geomFace, uv.X(), uv.Y() );
3539         }
3540         else
3541         {
3542           getMeshDS()->SetNodeInVolume( node, helper.GetSubShapeID() );
3543         }
3544       }
3545       else
3546       {
3547         if ( !edge._sWOL.IsNull() )
3548         {
3549           // make average pos from new and current parameters
3550           if ( isOnEdge )
3551           {
3552             u = 0.5 * ( u + helper.GetNodeU( geomEdge, node ));
3553             pos = curve->Value( u ).Transformed(loc);
3554           }
3555           else
3556           {
3557             uv = 0.5 * ( uv + helper.GetNodeUV( geomFace, node ));
3558             pos = surface->Value( uv.X(), uv.Y()).Transformed(loc);
3559           }
3560         }
3561         node->setXYZ( pos.X(), pos.Y(), pos.Z() );
3562       }
3563     }
3564   }
3565
3566   if ( !getMeshDS()->IsEmbeddedMode() )
3567     // Log node movement
3568     for ( size_t i = 0; i < data._edges.size(); ++i )
3569     {
3570       _LayerEdge& edge = *data._edges[i];
3571       SMESH_TNodeXYZ p ( edge._nodes.back() );
3572       getMeshDS()->MoveNode( p._node, p.X(), p.Y(), p.Z() );
3573     }
3574
3575   // TODO: make quadratic prisms and polyhedrons(?)
3576
3577   helper.SetElementsOnShape(true);
3578
3579   TopExp_Explorer exp( data._solid, TopAbs_FACE );
3580   for ( ; exp.More(); exp.Next() )
3581   {
3582     if ( data._ignoreFaceIds.count( getMeshDS()->ShapeToIndex( exp.Current() )))
3583       continue;
3584     SMESHDS_SubMesh* fSubM = getMeshDS()->MeshElements( exp.Current() );
3585     SMDS_ElemIteratorPtr fIt = fSubM->GetElements();
3586     vector< vector<const SMDS_MeshNode*>* > nnVec;
3587     while ( fIt->more() )
3588     {
3589       const SMDS_MeshElement* face = fIt->next();
3590       int nbNodes = face->NbCornerNodes();
3591       nnVec.resize( nbNodes );
3592       SMDS_ElemIteratorPtr nIt = face->nodesIterator();
3593       for ( int iN = 0; iN < nbNodes; ++iN )
3594       {
3595         const SMDS_MeshNode* n = static_cast<const SMDS_MeshNode*>( nIt->next() );
3596         nnVec[ iN ] = & data._n2eMap[ n ]->_nodes;
3597       }
3598
3599       int nbZ = nnVec[0]->size();
3600       switch ( nbNodes )
3601       {
3602       case 3:
3603         for ( int iZ = 1; iZ < nbZ; ++iZ )
3604           helper.AddVolume( (*nnVec[0])[iZ-1], (*nnVec[1])[iZ-1], (*nnVec[2])[iZ-1],
3605                             (*nnVec[0])[iZ],   (*nnVec[1])[iZ],   (*nnVec[2])[iZ]);
3606         break;
3607       case 4:
3608         for ( int iZ = 1; iZ < nbZ; ++iZ )
3609           helper.AddVolume( (*nnVec[0])[iZ-1], (*nnVec[1])[iZ-1],
3610                             (*nnVec[2])[iZ-1], (*nnVec[3])[iZ-1],
3611                             (*nnVec[0])[iZ],   (*nnVec[1])[iZ],
3612                             (*nnVec[2])[iZ],   (*nnVec[3])[iZ]);
3613         break;
3614       default:
3615         return error("Not supported type of element", data._index);
3616       }
3617     }
3618   }
3619   return true;
3620 }
3621
3622 //================================================================================
3623 /*!
3624  * \brief Shrink 2D mesh on faces to let space for inflated layers
3625  */
3626 //================================================================================
3627
3628 bool _ViscousBuilder::shrink()
3629 {
3630   // make map of (ids of FACEs to shrink mesh on) to (_SolidData containing _LayerEdge's
3631   // inflated along FACE or EDGE)
3632   map< TGeomID, _SolidData* > f2sdMap;
3633   for ( size_t i = 0 ; i < _sdVec.size(); ++i )
3634   {
3635     _SolidData& data = _sdVec[i];
3636     TopTools_MapOfShape FFMap;
3637     map< TGeomID, TopoDS_Shape >::iterator s2s = data._shrinkShape2Shape.begin();
3638     for (; s2s != data._shrinkShape2Shape.end(); ++s2s )
3639       if ( s2s->second.ShapeType() == TopAbs_FACE )
3640       {
3641         f2sdMap.insert( make_pair( getMeshDS()->ShapeToIndex( s2s->second ), &data ));
3642
3643         if ( FFMap.Add( (*s2s).second ))
3644           // Put mesh faces on the shrinked FACE to the proxy sub-mesh to avoid
3645           // usage of mesh faces made in addBoundaryElements() by the 3D algo or
3646           // by StdMeshers_QuadToTriaAdaptor
3647           if ( SMESHDS_SubMesh* smDS = getMeshDS()->MeshElements( s2s->second ))
3648           {
3649             SMESH_ProxyMesh::SubMesh* proxySub =
3650               data._proxyMesh->getFaceSubM( TopoDS::Face( s2s->second ), /*create=*/true);
3651             SMDS_ElemIteratorPtr fIt = smDS->GetElements();
3652             while ( fIt->more() )
3653               proxySub->AddElement( fIt->next() );
3654             // as a result 3D algo will use elements from proxySub and not from smDS
3655           }
3656       }
3657   }
3658
3659   SMESH_MesherHelper helper( *_mesh );
3660   helper.ToFixNodeParameters( true );
3661
3662   // EDGE's to shrink
3663   map< TGeomID, _Shrinker1D > e2shrMap;
3664
3665   // loop on FACES to srink mesh on
3666   map< TGeomID, _SolidData* >::iterator f2sd = f2sdMap.begin();
3667   for ( ; f2sd != f2sdMap.end(); ++f2sd )
3668   {
3669     _SolidData&     data = *f2sd->second;
3670     TNode2Edge&   n2eMap = data._n2eMap;
3671     const TopoDS_Face& F = TopoDS::Face( getMeshDS()->IndexToShape( f2sd->first ));
3672
3673     Handle(Geom_Surface) surface = BRep_Tool::Surface(F);
3674
3675     SMESH_subMesh*     sm = _mesh->GetSubMesh( F );
3676     SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
3677
3678     helper.SetSubShape(F);
3679
3680     // ===========================
3681     // Prepare data for shrinking
3682     // ===========================
3683
3684     // Collect nodes to smooth, as src nodes are not yet replaced by tgt ones
3685     // and thus all nodes on a FACE connected to 2d elements are to be smoothed
3686     vector < const SMDS_MeshNode* > smoothNodes;
3687     {
3688       SMDS_NodeIteratorPtr nIt = smDS->GetNodes();
3689       while ( nIt->more() )
3690       {
3691         const SMDS_MeshNode* n = nIt->next();
3692         if ( n->NbInverseElements( SMDSAbs_Face ) > 0 )
3693           smoothNodes.push_back( n );
3694       }
3695     }
3696     // Find out face orientation
3697     double refSign = 1;
3698     const set<TGeomID> ignoreShapes;
3699     bool isOkUV;
3700     if ( !smoothNodes.empty() )
3701     {
3702       vector<_Simplex> simplices;
3703       getSimplices( smoothNodes[0], simplices, ignoreShapes );
3704       helper.GetNodeUV( F, simplices[0]._nPrev, 0, &isOkUV ); // fix UV of silpmex nodes
3705       helper.GetNodeUV( F, simplices[0]._nNext, 0, &isOkUV );
3706       gp_XY uv = helper.GetNodeUV( F, smoothNodes[0], 0, &isOkUV );
3707       if ( !simplices[0].IsForward(uv, smoothNodes[0], F, helper,refSign) )
3708         refSign = -1;
3709     }
3710
3711     // Find _LayerEdge's inflated along F
3712     vector< _LayerEdge* > lEdges;
3713     {
3714       SMESH_subMeshIteratorPtr subIt =
3715         sm->getDependsOnIterator(/*includeSelf=*/false, /*complexShapeFirst=*/false);
3716       while ( subIt->more() )
3717       {
3718         SMESH_subMesh*     sub = subIt->next();
3719         SMESHDS_SubMesh* subDS = sub->GetSubMeshDS();
3720         if ( subDS->NbNodes() == 0 || !n2eMap.count( subDS->GetNodes()->next() ))
3721           continue;
3722         SMDS_NodeIteratorPtr nIt = subDS->GetNodes();
3723         while ( nIt->more() )
3724         {
3725           _LayerEdge* edge = n2eMap[ nIt->next() ];
3726           lEdges.push_back( edge );
3727           prepareEdgeToShrink( *edge, F, helper, smDS );
3728         }
3729       }
3730     }
3731
3732     dumpFunction(SMESH_Comment("beforeShrinkFace")<<f2sd->first); // debug
3733     SMDS_ElemIteratorPtr fIt = smDS->GetElements();
3734     while ( fIt->more() )
3735       if ( const SMDS_MeshElement* f = fIt->next() )
3736         dumpChangeNodes( f );
3737
3738     // Replace source nodes by target nodes in mesh faces to shrink
3739     const SMDS_MeshNode* nodes[20];
3740     for ( size_t i = 0; i < lEdges.size(); ++i )
3741     {
3742       _LayerEdge& edge = *lEdges[i];
3743       const SMDS_MeshNode* srcNode = edge._nodes[0];
3744       const SMDS_MeshNode* tgtNode = edge._nodes.back();
3745       SMDS_ElemIteratorPtr fIt = srcNode->GetInverseElementIterator(SMDSAbs_Face);
3746       while ( fIt->more() )
3747       {
3748         const SMDS_MeshElement* f = fIt->next();
3749         if ( !smDS->Contains( f ))
3750           continue;
3751         SMDS_NodeIteratorPtr nIt = f->nodeIterator();
3752         for ( int iN = 0; nIt->more(); ++iN )
3753         {
3754           const SMDS_MeshNode* n = nIt->next();
3755           nodes[iN] = ( n == srcNode ? tgtNode : n );
3756         }
3757         helper.GetMeshDS()->ChangeElementNodes( f, nodes, f->NbNodes() );
3758       }
3759     }
3760
3761     // find out if a FACE is concave
3762     const bool isConcaveFace = isConcave( F, helper );
3763
3764     // Create _SmoothNode's on face F
3765     vector< _SmoothNode > nodesToSmooth( smoothNodes.size() );
3766     {
3767       const bool sortSimplices = isConcaveFace;
3768       for ( size_t i = 0; i < smoothNodes.size(); ++i )
3769       {
3770         const SMDS_MeshNode* n = smoothNodes[i];
3771         nodesToSmooth[ i ]._node = n;
3772         // src nodes must be replaced by tgt nodes to have tgt nodes in _simplices
3773         getSimplices( n, nodesToSmooth[ i ]._simplices, ignoreShapes, NULL, sortSimplices );
3774         // fix up incorrect uv of nodes on the FACE
3775         helper.GetNodeUV( F, n, 0, &isOkUV);
3776         dumpMove( n );
3777       }
3778     }
3779     //if ( nodesToSmooth.empty() ) continue;
3780
3781     // Find EDGE's to shrink and set simpices to LayerEdge's
3782     set< _Shrinker1D* > eShri1D;
3783     {
3784       for ( size_t i = 0; i < lEdges.size(); ++i )
3785       {
3786         _LayerEdge* edge = lEdges[i];
3787         if ( edge->_sWOL.ShapeType() == TopAbs_EDGE )
3788         {
3789           TGeomID edgeIndex = getMeshDS()->ShapeToIndex( edge->_sWOL );
3790           _Shrinker1D& srinker = e2shrMap[ edgeIndex ];
3791           eShri1D.insert( & srinker );
3792           srinker.AddEdge( edge, helper );
3793           VISCOUS_3D::ToClearSubWithMain( _mesh->GetSubMesh( edge->_sWOL ), data._solid );
3794           // restore params of nodes on EGDE if the EDGE has been already
3795           // srinked while srinking another FACE
3796           srinker.RestoreParams();
3797         }
3798         getSimplices( /*tgtNode=*/edge->_nodes.back(), edge->_simplices, ignoreShapes );
3799       }
3800     }
3801
3802     bool toFixTria = false; // to improve quality of trias by diagonal swap
3803     if ( isConcaveFace )
3804     {
3805       const bool hasTria = _mesh->NbTriangles(), hasQuad = _mesh->NbQuadrangles();
3806       if ( hasTria != hasQuad ) {
3807         toFixTria = hasTria;
3808       }
3809       else {
3810         set<int> nbNodesSet;
3811         SMDS_ElemIteratorPtr fIt = smDS->GetElements();
3812         while ( fIt->more() && nbNodesSet.size() < 2 )
3813           nbNodesSet.insert( fIt->next()->NbCornerNodes() );
3814         toFixTria = ( *nbNodesSet.begin() == 3 );
3815       }
3816     }
3817
3818     // ==================
3819     // Perform shrinking
3820     // ==================
3821
3822     bool shrinked = true;
3823     int badNb, shriStep=0, smooStep=0;
3824     _SmoothNode::SmoothType smoothType
3825       = isConcaveFace ? _SmoothNode::ANGULAR : _SmoothNode::LAPLACIAN;
3826     while ( shrinked )
3827     {
3828       shriStep++;
3829       // Move boundary nodes (actually just set new UV)
3830       // -----------------------------------------------
3831       dumpFunction(SMESH_Comment("moveBoundaryOnF")<<f2sd->first<<"_st"<<shriStep ); // debug
3832       shrinked = false;
3833       for ( size_t i = 0; i < lEdges.size(); ++i )
3834       {
3835         shrinked |= lEdges[i]->SetNewLength2d( surface,F,helper );
3836       }
3837       dumpFunctionEnd();
3838
3839       // Move nodes on EDGE's
3840       // (XYZ is set as soon as a needed length reached in SetNewLength2d())
3841       set< _Shrinker1D* >::iterator shr = eShri1D.begin();
3842       for ( ; shr != eShri1D.end(); ++shr )
3843         (*shr)->Compute( /*set3D=*/false, helper );
3844
3845       // Smoothing in 2D
3846       // -----------------
3847       int nbNoImpSteps = 0;
3848       bool       moved = true;
3849       badNb = 1;
3850       while (( nbNoImpSteps < 5 && badNb > 0) && moved)
3851       {
3852         dumpFunction(SMESH_Comment("shrinkFace")<<f2sd->first<<"_st"<<++smooStep); // debug
3853
3854         int oldBadNb = badNb;
3855         badNb = 0;
3856         moved = false;
3857         for ( size_t i = 0; i < nodesToSmooth.size(); ++i )
3858         {
3859           moved |= nodesToSmooth[i].Smooth( badNb,surface,helper,refSign,
3860                                             smoothType, /*set3D=*/isConcaveFace);
3861         }
3862         if ( badNb < oldBadNb )
3863           nbNoImpSteps = 0;
3864         else
3865           nbNoImpSteps++;
3866
3867         dumpFunctionEnd();
3868       }
3869       if ( badNb > 0 )
3870         return error(SMESH_Comment("Can't shrink 2D mesh on face ") << f2sd->first );
3871       if ( shriStep > 200 )
3872         return error(SMESH_Comment("Infinite loop at shrinking 2D mesh on face ") << f2sd->first );
3873
3874       // Fix narrow triangles by swapping diagonals
3875       // ---------------------------------------
3876       if ( toFixTria )
3877       {
3878         set<const SMDS_MeshNode*> usedNodes;
3879         fixBadFaces( F, helper, /*is2D=*/true, shriStep, & usedNodes); // swap diagonals
3880
3881         // update working data
3882         set<const SMDS_MeshNode*>::iterator n;
3883         for ( size_t i = 0; i < nodesToSmooth.size() && !usedNodes.empty(); ++i )
3884         {
3885           n = usedNodes.find( nodesToSmooth[ i ]._node );
3886           if ( n != usedNodes.end())
3887           {
3888             getSimplices( nodesToSmooth[ i ]._node,
3889                           nodesToSmooth[ i ]._simplices,
3890                           ignoreShapes, NULL,
3891                           /*sortSimplices=*/ smoothType == _SmoothNode::ANGULAR );
3892             usedNodes.erase( n );
3893           }
3894         }
3895         for ( size_t i = 0; i < lEdges.size() && !usedNodes.empty(); ++i )
3896         {
3897           n = usedNodes.find( /*tgtNode=*/ lEdges[i]->_nodes.back() );
3898           if ( n != usedNodes.end())
3899           {
3900             getSimplices( lEdges[i]->_nodes.back(),
3901                           lEdges[i]->_simplices,
3902                           ignoreShapes );
3903             usedNodes.erase( n );
3904           }
3905         }
3906       }
3907     } // while ( shrinked )
3908
3909     // No wrongly shaped faces remain; final smooth. Set node XYZ.
3910     bool isStructuredFixed = false;
3911     if ( SMESH_2D_Algo* algo = dynamic_cast<SMESH_2D_Algo*>( sm->GetAlgo() ))
3912       isStructuredFixed = algo->FixInternalNodes( *data._proxyMesh, F );
3913     if ( !isStructuredFixed )
3914     {
3915       if ( isConcaveFace ) // fix narrow faces by swapping diagonals
3916         fixBadFaces( F, helper, /*is2D=*/false, ++shriStep );
3917
3918       for ( int st = 3; st; --st )
3919       {
3920         switch( st ) {
3921         case 1: smoothType = _SmoothNode::LAPLACIAN; break;
3922         case 2: smoothType = _SmoothNode::LAPLACIAN; break;
3923         case 3: smoothType = _SmoothNode::ANGULAR; break;
3924         }
3925         dumpFunction(SMESH_Comment("shrinkFace")<<f2sd->first<<"_st"<<++smooStep); // debug
3926         for ( size_t i = 0; i < nodesToSmooth.size(); ++i )
3927         {
3928           nodesToSmooth[i].Smooth( badNb,surface,helper,refSign,
3929                                    smoothType,/*set3D=*/st==1 );
3930         }
3931         dumpFunctionEnd();
3932       }
3933     }
3934     // Set an event listener to clear FACE sub-mesh together with SOLID sub-mesh
3935     VISCOUS_3D::ToClearSubWithMain( sm, data._solid );
3936
3937     if ( !getMeshDS()->IsEmbeddedMode() )
3938       // Log node movement
3939       for ( size_t i = 0; i < nodesToSmooth.size(); ++i )
3940       {
3941         SMESH_TNodeXYZ p ( nodesToSmooth[i]._node );
3942         getMeshDS()->MoveNode( nodesToSmooth[i]._node, p.X(), p.Y(), p.Z() );
3943       }
3944
3945   } // loop on FACES to srink mesh on
3946
3947
3948   // Replace source nodes by target nodes in shrinked mesh edges
3949
3950   map< int, _Shrinker1D >::iterator e2shr = e2shrMap.begin();
3951   for ( ; e2shr != e2shrMap.end(); ++e2shr )
3952     e2shr->second.SwapSrcTgtNodes( getMeshDS() );
3953
3954   return true;
3955 }
3956
3957 //================================================================================
3958 /*!
3959  * \brief Computes 2d shrink direction and finds nodes limiting shrinking
3960  */
3961 //================================================================================
3962
3963 bool _ViscousBuilder::prepareEdgeToShrink( _LayerEdge&            edge,
3964                                            const TopoDS_Face&     F,
3965                                            SMESH_MesherHelper&    helper,
3966                                            const SMESHDS_SubMesh* faceSubMesh)
3967 {
3968   const SMDS_MeshNode* srcNode = edge._nodes[0];
3969   const SMDS_MeshNode* tgtNode = edge._nodes.back();
3970
3971   edge._pos.clear();
3972
3973   if ( edge._sWOL.ShapeType() == TopAbs_FACE )
3974   {
3975     gp_XY srcUV = helper.GetNodeUV( F, srcNode );
3976     gp_XY tgtUV = helper.GetNodeUV( F, tgtNode );
3977     gp_Vec2d uvDir( srcUV, tgtUV );
3978     double uvLen = uvDir.Magnitude();
3979     uvDir /= uvLen;
3980     edge._normal.SetCoord( uvDir.X(),uvDir.Y(), 0);
3981     edge._len = uvLen;
3982
3983     // // IMPORTANT to have src nodes NOT yet REPLACED by tgt nodes in shrinked faces
3984     // vector<const SMDS_MeshElement*> faces;
3985     // multimap< double, const SMDS_MeshNode* > proj2node;
3986     // SMDS_ElemIteratorPtr fIt = srcNode->GetInverseElementIterator(SMDSAbs_Face);
3987     // while ( fIt->more() )
3988     // {
3989     //   const SMDS_MeshElement* f = fIt->next();
3990     //   if ( faceSubMesh->Contains( f ))
3991     //     faces.push_back( f );
3992     // }
3993     // for ( size_t i = 0; i < faces.size(); ++i )
3994     // {
3995     //   const int nbNodes = faces[i]->NbCornerNodes();
3996     //   for ( int j = 0; j < nbNodes; ++j )
3997     //   {
3998     //     const SMDS_MeshNode* n = faces[i]->GetNode(j);
3999     //     if ( n == srcNode ) continue;
4000     //     if ( n->GetPosition()->GetTypeOfPosition() != SMDS_TOP_FACE &&
4001     //          ( faces.size() > 1 || nbNodes > 3 ))
4002     //       continue;
4003     //     gp_Pnt2d uv = helper.GetNodeUV( F, n );
4004     //     gp_Vec2d uvDirN( srcUV, uv );
4005     //     double proj = uvDirN * uvDir;
4006     //     proj2node.insert( make_pair( proj, n ));
4007     //   }
4008     // }
4009
4010     // multimap< double, const SMDS_MeshNode* >::iterator p2n = proj2node.begin(), p2nEnd;
4011     // const double       minProj = p2n->first;
4012     // const double projThreshold = 1.1 * uvLen;
4013     // if ( minProj > projThreshold )
4014     // {
4015     //   // tgtNode is located so that it does not make faces with wrong orientation
4016     //   return true;
4017     // }
4018     edge._pos.resize(1);
4019     edge._pos[0].SetCoord( tgtUV.X(), tgtUV.Y(), 0 );
4020
4021     // store most risky nodes in _simplices
4022     // p2nEnd = proj2node.lower_bound( projThreshold );
4023     // int nbSimpl = ( std::distance( p2n, p2nEnd ) + 1) / 2;
4024     // edge._simplices.resize( nbSimpl );
4025     // for ( int i = 0; i < nbSimpl; ++i )
4026     // {
4027     //   edge._simplices[i]._nPrev = p2n->second;
4028     //   if ( ++p2n != p2nEnd )
4029     //     edge._simplices[i]._nNext = p2n->second;
4030     // }
4031     // set UV of source node to target node
4032     SMDS_FacePosition* pos = static_cast<SMDS_FacePosition*>( tgtNode->GetPosition() );
4033     pos->SetUParameter( srcUV.X() );
4034     pos->SetVParameter( srcUV.Y() );
4035   }
4036   else // _sWOL is TopAbs_EDGE
4037   {
4038     TopoDS_Edge E = TopoDS::Edge( edge._sWOL);
4039     SMESHDS_SubMesh* edgeSM = getMeshDS()->MeshElements( E );
4040     if ( !edgeSM || edgeSM->NbElements() == 0 )
4041       return error(SMESH_Comment("Not meshed EDGE ") << getMeshDS()->ShapeToIndex( E ));
4042
4043     const SMDS_MeshNode* n2 = 0;
4044     SMDS_ElemIteratorPtr eIt = srcNode->GetInverseElementIterator(SMDSAbs_Edge);
4045     while ( eIt->more() && !n2 )
4046     {
4047       const SMDS_MeshElement* e = eIt->next();
4048       if ( !edgeSM->Contains(e)) continue;
4049       n2 = e->GetNode( 0 );
4050       if ( n2 == srcNode ) n2 = e->GetNode( 1 );
4051     }
4052     if ( !n2 )
4053       return error(SMESH_Comment("Wrongly meshed EDGE ") << getMeshDS()->ShapeToIndex( E ));
4054
4055     double uSrc = helper.GetNodeU( E, srcNode, n2 );
4056     double uTgt = helper.GetNodeU( E, tgtNode, srcNode );
4057     double u2   = helper.GetNodeU( E, n2,      srcNode );
4058
4059     if ( fabs( uSrc-uTgt ) < 0.99 * fabs( uSrc-u2 ))
4060     {
4061       // tgtNode is located so that it does not make faces with wrong orientation
4062       return true;
4063     }
4064     edge._pos.resize(1);
4065     edge._pos[0].SetCoord( U_TGT, uTgt );
4066     edge._pos[0].SetCoord( U_SRC, uSrc );
4067     edge._pos[0].SetCoord( LEN_TGT, fabs( uSrc-uTgt ));
4068
4069     edge._simplices.resize( 1 );
4070     edge._simplices[0]._nPrev = n2;
4071
4072     // set UV of source node to target node
4073     SMDS_EdgePosition* pos = static_cast<SMDS_EdgePosition*>( tgtNode->GetPosition() );
4074     pos->SetUParameter( uSrc );
4075   }
4076   return true;
4077
4078   //================================================================================
4079   /*!
4080    * \brief Compute positions (UV) to set to a node on edge moved during shrinking
4081    */
4082   //================================================================================
4083   
4084   // Compute UV to follow during shrinking
4085
4086 //   const SMDS_MeshNode* srcNode = edge._nodes[0];
4087 //   const SMDS_MeshNode* tgtNode = edge._nodes.back();
4088
4089 //   gp_XY srcUV = helper.GetNodeUV( F, srcNode );
4090 //   gp_XY tgtUV = helper.GetNodeUV( F, tgtNode );
4091 //   gp_Vec2d uvDir( srcUV, tgtUV );
4092 //   double uvLen = uvDir.Magnitude();
4093 //   uvDir /= uvLen;
4094
4095 //   // Select shrinking step such that not to make faces with wrong orientation.
4096 //   // IMPORTANT to have src nodes NOT yet REPLACED by tgt nodes in shrinked faces
4097 //   const double minStepSize = uvLen / 20;
4098 //   double stepSize = uvLen;
4099 //   SMDS_ElemIteratorPtr fIt = srcNode->GetInverseElementIterator(SMDSAbs_Face);
4100 //   while ( fIt->more() )
4101 //   {
4102 //     const SMDS_MeshElement* f = fIt->next();
4103 //     if ( !faceSubMesh->Contains( f )) continue;
4104 //     const int nbNodes = f->NbCornerNodes();
4105 //     for ( int i = 0; i < nbNodes; ++i )
4106 //     {
4107 //       const SMDS_MeshNode* n = f->GetNode(i);
4108 //       if ( n->GetPosition()->GetTypeOfPosition() != SMDS_TOP_FACE || n == srcNode)
4109 //         continue;
4110 //       gp_XY uv = helper.GetNodeUV( F, n );
4111 //       gp_Vec2d uvDirN( srcUV, uv );
4112 //       double proj = uvDirN * uvDir;
4113 //       if ( proj < stepSize && proj > minStepSize )
4114 //         stepSize = proj;
4115 //     }
4116 //   }
4117 //   stepSize *= 0.8;
4118
4119 //   const int nbSteps = ceil( uvLen / stepSize );
4120 //   gp_XYZ srcUV0( srcUV.X(), srcUV.Y(), 0 );
4121 //   gp_XYZ tgtUV0( tgtUV.X(), tgtUV.Y(), 0 );
4122 //   edge._pos.resize( nbSteps );
4123 //   edge._pos[0] = tgtUV0;
4124 //   for ( int i = 1; i < nbSteps; ++i )
4125 //   {
4126 //     double r = i / double( nbSteps );
4127 //     edge._pos[i] = (1-r) * tgtUV0 + r * srcUV0;
4128 //   }
4129 //   return true;
4130 }
4131
4132 //================================================================================
4133 /*!
4134  * \brief Try to fix triangles with high aspect ratio by swaping diagonals
4135  */
4136 //================================================================================
4137
4138 void _ViscousBuilder::fixBadFaces(const TopoDS_Face&          F,
4139                                   SMESH_MesherHelper&         helper,
4140                                   const bool                  is2D,
4141                                   const int                   step,
4142                                   set<const SMDS_MeshNode*> * involvedNodes)
4143 {
4144   SMESH::Controls::AspectRatio qualifier;
4145   SMESH::Controls::TSequenceOfXYZ points(3), points1(3), points2(3);
4146   const double maxAspectRatio = is2D ? 4. : 2;
4147   NodeCoordHelper xyz( F, helper, is2D );
4148
4149   // find bad triangles
4150
4151   vector< const SMDS_MeshElement* > badTrias;
4152   vector< double >                  badAspects;
4153   SMESHDS_SubMesh*      sm = helper.GetMeshDS()->MeshElements( F );
4154   SMDS_ElemIteratorPtr fIt = sm->GetElements();
4155   while ( fIt->more() )
4156   {
4157     const SMDS_MeshElement * f = fIt->next();
4158     if ( f->NbCornerNodes() != 3 ) continue;
4159     for ( int iP = 0; iP < 3; ++iP ) points(iP+1) = xyz( f->GetNode(iP));
4160     double aspect = qualifier.GetValue( points );
4161     if ( aspect > maxAspectRatio )
4162     {
4163       badTrias.push_back( f );
4164       badAspects.push_back( aspect );
4165     }
4166   }
4167   if ( step == 1 )
4168   {
4169     dumpFunction(SMESH_Comment("beforeSwapDiagonals_F")<<helper.GetSubShapeID());
4170     SMDS_ElemIteratorPtr fIt = sm->GetElements();
4171     while ( fIt->more() )
4172     {
4173       const SMDS_MeshElement * f = fIt->next();
4174       if ( f->NbCornerNodes() == 3 )
4175         dumpChangeNodes( f );
4176     }
4177     dumpFunctionEnd();
4178   }
4179   if ( badTrias.empty() )
4180     return;
4181
4182   // find couples of faces to swap diagonal
4183
4184   typedef pair < const SMDS_MeshElement* , const SMDS_MeshElement* > T2Trias;
4185   vector< T2Trias > triaCouples; 
4186
4187   TIDSortedElemSet involvedFaces, emptySet;
4188   for ( size_t iTia = 0; iTia < badTrias.size(); ++iTia )
4189   {
4190     T2Trias trias    [3];
4191     double  aspRatio [3];
4192     int i1, i2, i3;
4193
4194     if ( !involvedFaces.insert( badTrias[iTia] ).second )
4195       continue;
4196     for ( int iP = 0; iP < 3; ++iP )
4197       points(iP+1) = xyz( badTrias[iTia]->GetNode(iP));
4198
4199     // find triangles adjacent to badTrias[iTia] with better aspect ratio after diag-swaping
4200     int bestCouple = -1;
4201     for ( int iSide = 0; iSide < 3; ++iSide )
4202     {
4203       const SMDS_MeshNode* n1 = badTrias[iTia]->GetNode( iSide );
4204       const SMDS_MeshNode* n2 = badTrias[iTia]->GetNode(( iSide+1 ) % 3 );
4205       trias [iSide].first  = badTrias[iTia];
4206       trias [iSide].second = SMESH_MeshAlgos::FindFaceInSet( n1, n2, emptySet, involvedFaces,
4207                                                              & i1, & i2 );
4208       if ( ! trias[iSide].second || trias[iSide].second->NbCornerNodes() != 3 )
4209         continue;
4210
4211       // aspect ratio of an adjacent tria
4212       for ( int iP = 0; iP < 3; ++iP )
4213         points2(iP+1) = xyz( trias[iSide].second->GetNode(iP));
4214       double aspectInit = qualifier.GetValue( points2 );
4215
4216       // arrange nodes as after diag-swaping
4217       if ( helper.WrapIndex( i1+1, 3 ) == i2 )
4218         i3 = helper.WrapIndex( i1-1, 3 );
4219       else
4220         i3 = helper.WrapIndex( i1+1, 3 );
4221       points1 = points;
4222       points1( 1+ iSide ) = points2( 1+ i3 );
4223       points2( 1+ i2    ) = points1( 1+ ( iSide+2 ) % 3 );
4224
4225       // aspect ratio after diag-swaping
4226       aspRatio[ iSide ] = qualifier.GetValue( points1 ) + qualifier.GetValue( points2 );
4227       if ( aspRatio[ iSide ] > aspectInit + badAspects[ iTia ] )
4228         continue;
4229
4230       // prevent inversion of a triangle
4231       gp_Vec norm1 = gp_Vec( points1(1), points1(3) ) ^ gp_Vec( points1(1), points1(2) );
4232       gp_Vec norm2 = gp_Vec( points2(1), points2(3) ) ^ gp_Vec( points2(1), points2(2) );
4233       if ( norm1 * norm2 < 0. && norm1.Angle( norm2 ) > 70./180.*M_PI )
4234         continue;
4235
4236       if ( bestCouple < 0 || aspRatio[ bestCouple ] > aspRatio[ iSide ] )
4237         bestCouple = iSide;
4238     }
4239
4240     if ( bestCouple >= 0 )
4241     {
4242       triaCouples.push_back( trias[bestCouple] );
4243       involvedFaces.insert ( trias[bestCouple].second );
4244     }
4245     else
4246     {
4247       involvedFaces.erase( badTrias[iTia] );
4248     }
4249   }
4250   if ( triaCouples.empty() )
4251     return;
4252
4253   // swap diagonals
4254
4255   SMESH_MeshEditor editor( helper.GetMesh() );
4256   dumpFunction(SMESH_Comment("beforeSwapDiagonals_F")<<helper.GetSubShapeID()<<"_"<<step);
4257   for ( size_t i = 0; i < triaCouples.size(); ++i )
4258   {
4259     dumpChangeNodes( triaCouples[i].first );
4260     dumpChangeNodes( triaCouples[i].second );
4261     editor.InverseDiag( triaCouples[i].first, triaCouples[i].second );
4262   }
4263
4264   if ( involvedNodes )
4265     for ( size_t i = 0; i < triaCouples.size(); ++i )
4266     {
4267       involvedNodes->insert( triaCouples[i].first->begin_nodes(),
4268                              triaCouples[i].first->end_nodes() );
4269       involvedNodes->insert( triaCouples[i].second->begin_nodes(),
4270                              triaCouples[i].second->end_nodes() );
4271     }
4272
4273   // just for debug dump resulting triangles
4274   dumpFunction(SMESH_Comment("swapDiagonals_F")<<helper.GetSubShapeID()<<"_"<<step);
4275   for ( size_t i = 0; i < triaCouples.size(); ++i )
4276   {
4277     dumpChangeNodes( triaCouples[i].first );
4278     dumpChangeNodes( triaCouples[i].second );
4279   }
4280 }
4281
4282 //================================================================================
4283 /*!
4284  * \brief Move target node to it's final position on the FACE during shrinking
4285  */
4286 //================================================================================
4287
4288 bool _LayerEdge::SetNewLength2d( Handle(Geom_Surface)& surface,
4289                                  const TopoDS_Face&    F,
4290                                  SMESH_MesherHelper&   helper )
4291 {
4292   if ( _pos.empty() )
4293     return false; // already at the target position
4294
4295   SMDS_MeshNode* tgtNode = const_cast< SMDS_MeshNode*& >( _nodes.back() );
4296
4297   if ( _sWOL.ShapeType() == TopAbs_FACE )
4298   {
4299     gp_XY    curUV = helper.GetNodeUV( F, tgtNode );
4300     gp_Pnt2d tgtUV( _pos[0].X(), _pos[0].Y() );
4301     gp_Vec2d uvDir( _normal.X(), _normal.Y() );
4302     const double uvLen = tgtUV.Distance( curUV );
4303     const double kSafe = Max( 0.5, 1. - 0.1 * _simplices.size() );
4304
4305     // Select shrinking step such that not to make faces with wrong orientation.
4306     double stepSize = uvLen;
4307     for ( size_t i = 0; i < _simplices.size(); ++i )
4308     {
4309       // find intersection of 2 lines: curUV-tgtUV and that connecting simplex nodes
4310       gp_XY uvN1 = helper.GetNodeUV( F, _simplices[i]._nPrev );
4311       gp_XY uvN2 = helper.GetNodeUV( F, _simplices[i]._nNext );
4312       gp_XY dirN = uvN2 - uvN1;
4313       double det = uvDir.Crossed( dirN );
4314       if ( Abs( det )  < std::numeric_limits<double>::min() ) continue;
4315       gp_XY dirN2Cur = curUV - uvN1;
4316       double step = dirN.Crossed( dirN2Cur ) / det;
4317       if ( step > 0 )
4318         stepSize = Min( step, stepSize );
4319     }
4320     gp_Pnt2d newUV;
4321     if ( uvLen - stepSize < _len / 200. )
4322     {
4323       newUV = tgtUV;
4324       _pos.clear();
4325     }
4326     else if ( stepSize > 0 )
4327     {
4328       newUV = curUV + uvDir.XY() * stepSize * kSafe;
4329     }
4330     else
4331     {
4332       return true;
4333     }
4334     SMDS_FacePosition* pos = static_cast<SMDS_FacePosition*>( tgtNode->GetPosition() );
4335     pos->SetUParameter( newUV.X() );
4336     pos->SetVParameter( newUV.Y() );
4337
4338 #ifdef __myDEBUG
4339     gp_Pnt p = surface->Value( newUV.X(), newUV.Y() );
4340     tgtNode->setXYZ( p.X(), p.Y(), p.Z() );
4341     dumpMove( tgtNode );
4342 #endif
4343   }
4344   else // _sWOL is TopAbs_EDGE
4345   {
4346     TopoDS_Edge E = TopoDS::Edge( _sWOL );
4347     const SMDS_MeshNode* n2 = _simplices[0]._nPrev;
4348     SMDS_EdgePosition* tgtPos = static_cast<SMDS_EdgePosition*>( tgtNode->GetPosition() );
4349
4350     const double u2 = helper.GetNodeU( E, n2, tgtNode );
4351     const double uSrc   = _pos[0].Coord( U_SRC );
4352     const double lenTgt = _pos[0].Coord( LEN_TGT );
4353
4354     double newU = _pos[0].Coord( U_TGT );
4355     if ( lenTgt < 0.99 * fabs( uSrc-u2 )) // n2 got out of src-tgt range
4356     {
4357       _pos.clear();
4358     }
4359     else
4360     {
4361       newU = 0.1 * tgtPos->GetUParameter() + 0.9 * u2;
4362     }
4363     tgtPos->SetUParameter( newU );
4364 #ifdef __myDEBUG
4365     gp_XY newUV = helper.GetNodeUV( F, tgtNode, _nodes[0]);
4366     gp_Pnt p = surface->Value( newUV.X(), newUV.Y() );
4367     tgtNode->setXYZ( p.X(), p.Y(), p.Z() );
4368     dumpMove( tgtNode );
4369 #endif
4370   }
4371   return true;
4372 }
4373
4374 //================================================================================
4375 /*!
4376  * \brief Perform smooth on the FACE
4377  *  \retval bool - true if the node has been moved
4378  */
4379 //================================================================================
4380
4381 bool _SmoothNode::Smooth(int&                  badNb,
4382                          Handle(Geom_Surface)& surface,
4383                          SMESH_MesherHelper&   helper,
4384                          const double          refSign,
4385                          SmoothType            how,
4386                          bool                  set3D)
4387 {
4388   const TopoDS_Face& face = TopoDS::Face( helper.GetSubShape() );
4389
4390   // get uv of surrounding nodes
4391   vector<gp_XY> uv( _simplices.size() );
4392   for ( size_t i = 0; i < _simplices.size(); ++i )
4393     uv[i] = helper.GetNodeUV( face, _simplices[i]._nPrev, _node );
4394
4395   // compute new UV for the node
4396   gp_XY newPos (0,0);
4397   if ( how == TFI && _simplices.size() == 4 )
4398   {
4399     gp_XY corners[4];
4400     for ( size_t i = 0; i < _simplices.size(); ++i )
4401       if ( _simplices[i]._nOpp )
4402         corners[i] = helper.GetNodeUV( face, _simplices[i]._nOpp, _node );
4403       else
4404         throw SALOME_Exception(LOCALIZED("TFI smoothing: _Simplex::_nOpp not set!"));
4405
4406     newPos = helper.calcTFI ( 0.5, 0.5,
4407                               corners[0], corners[1], corners[2], corners[3],
4408                               uv[1], uv[2], uv[3], uv[0] );
4409   }
4410   else if ( how == ANGULAR )
4411   {
4412     newPos = computeAngularPos( uv, helper.GetNodeUV( face, _node ), refSign );
4413   }
4414   else if ( how == CENTROIDAL && _simplices.size() > 3 )
4415   {
4416     // average centers of diagonals wieghted with their reciprocal lengths
4417     if ( _simplices.size() == 4 )
4418     {
4419       double w1 = 1. / ( uv[2]-uv[0] ).SquareModulus();
4420       double w2 = 1. / ( uv[3]-uv[1] ).SquareModulus();
4421       newPos = ( w1 * ( uv[2]+uv[0] ) + w2 * ( uv[3]+uv[1] )) / ( w1+w2 ) / 2;
4422     }
4423     else
4424     {
4425       double sumWeight = 0;
4426       int nb = _simplices.size() == 4 ? 2 : _simplices.size();
4427       for ( int i = 0; i < nb; ++i )
4428       {
4429         int iFrom = i + 2;
4430         int iTo   = i + _simplices.size() - 1;
4431         for ( int j = iFrom; j < iTo; ++j )
4432         {
4433           int i2 = SMESH_MesherHelper::WrapIndex( j, _simplices.size() );
4434           double w = 1. / ( uv[i]-uv[i2] ).SquareModulus();
4435           sumWeight += w;
4436           newPos += w * ( uv[i]+uv[i2] );
4437         }
4438       }
4439       newPos /= 2 * sumWeight; // 2 is to get a middle between uv's
4440     }
4441   }
4442   else
4443   {
4444     // Laplacian smooth
4445     for ( size_t i = 0; i < _simplices.size(); ++i )
4446       newPos += uv[i];
4447     newPos /= _simplices.size();
4448   }
4449
4450   // count quality metrics (orientation) of triangles around the node
4451   int nbOkBefore = 0;
4452   gp_XY tgtUV = helper.GetNodeUV( face, _node );
4453   for ( size_t i = 0; i < _simplices.size(); ++i )
4454     nbOkBefore += _simplices[i].IsForward( tgtUV, _node, face, helper, refSign );
4455
4456   int nbOkAfter = 0;
4457   for ( size_t i = 0; i < _simplices.size(); ++i )
4458     nbOkAfter += _simplices[i].IsForward( newPos, _node, face, helper, refSign );
4459
4460   if ( nbOkAfter < nbOkBefore )
4461   {
4462     badNb += _simplices.size() - nbOkBefore;
4463     return false;
4464   }
4465
4466   SMDS_FacePosition* pos = static_cast<SMDS_FacePosition*>( _node->GetPosition() );
4467   pos->SetUParameter( newPos.X() );
4468   pos->SetVParameter( newPos.Y() );
4469
4470 #ifdef __myDEBUG
4471   set3D = true;
4472 #endif
4473   if ( set3D )
4474   {
4475     gp_Pnt p = surface->Value( newPos.X(), newPos.Y() );
4476     const_cast< SMDS_MeshNode* >( _node )->setXYZ( p.X(), p.Y(), p.Z() );
4477     dumpMove( _node );
4478   }
4479
4480   badNb += _simplices.size() - nbOkAfter;
4481   return ( (tgtUV-newPos).SquareModulus() > 1e-10 );
4482 }
4483
4484 //================================================================================
4485 /*!
4486  * \brief Computes new UV using angle based smoothing technic
4487  */
4488 //================================================================================
4489
4490 gp_XY _SmoothNode::computeAngularPos(vector<gp_XY>& uv,
4491                                      const gp_XY&   uvToFix,
4492                                      const double   refSign)
4493 {
4494   uv.push_back( uv.front() );
4495
4496   vector< gp_XY >  edgeDir ( uv.size() );
4497   vector< double > edgeSize( uv.size() );
4498   for ( size_t i = 1; i < edgeDir.size(); ++i )
4499   {
4500     edgeDir [i-1] = uv[i] - uv[i-1];
4501     edgeSize[i-1] = edgeDir[i-1].Modulus();
4502     if ( edgeSize[i-1] < numeric_limits<double>::min() )
4503       edgeDir[i-1].SetX( 100 );
4504     else
4505       edgeDir[i-1] /= edgeSize[i-1] * refSign;
4506   }
4507   edgeDir.back()  = edgeDir.front();
4508   edgeSize.back() = edgeSize.front();
4509
4510   gp_XY  newPos(0,0);
4511   int    nbEdges = 0;
4512   double sumSize = 0;
4513   for ( size_t i = 1; i < edgeDir.size(); ++i )
4514   {
4515     if ( edgeDir[i-1].X() > 1. ) continue;
4516     int i1 = i-1;
4517     while ( edgeDir[i].X() > 1. && ++i < edgeDir.size() );
4518     if ( i == edgeDir.size() ) break;
4519     gp_XY p = uv[i];
4520     gp_XY norm1( -edgeDir[i1].Y(), edgeDir[i1].X() );
4521     gp_XY norm2( -edgeDir[i].Y(),  edgeDir[i].X() );
4522     gp_XY bisec = norm1 + norm2;
4523     double bisecSize = bisec.Modulus();
4524     if ( bisecSize < numeric_limits<double>::min() )
4525     {
4526       bisec = -edgeDir[i1] + edgeDir[i];
4527       bisecSize = bisec.Modulus();
4528     }
4529     bisec /= bisecSize;
4530
4531     gp_XY  dirToN  = uvToFix - p;
4532     double distToN = dirToN.Modulus();
4533     if ( bisec * dirToN < 0 )
4534       distToN = -distToN;
4535
4536     newPos += ( p + bisec * distToN ) * ( edgeSize[i1] + edgeSize[i] );
4537     ++nbEdges;
4538     sumSize += edgeSize[i1] + edgeSize[i];
4539   }
4540   newPos /= /*nbEdges * */sumSize;
4541   return newPos;
4542 }
4543
4544 //================================================================================
4545 /*!
4546  * \brief Delete _SolidData
4547  */
4548 //================================================================================
4549
4550 _SolidData::~_SolidData()
4551 {
4552   for ( size_t i = 0; i < _edges.size(); ++i )
4553   {
4554     if ( _edges[i] && _edges[i]->_2neibors )
4555       delete _edges[i]->_2neibors;
4556     delete _edges[i];
4557   }
4558   _edges.clear();
4559 }
4560 //================================================================================
4561 /*!
4562  * \brief Add a _LayerEdge inflated along the EDGE
4563  */
4564 //================================================================================
4565
4566 void _Shrinker1D::AddEdge( const _LayerEdge* e, SMESH_MesherHelper& helper )
4567 {
4568   // init
4569   if ( _nodes.empty() )
4570   {
4571     _edges[0] = _edges[1] = 0;
4572     _done = false;
4573   }
4574   // check _LayerEdge
4575   if ( e == _edges[0] || e == _edges[1] )
4576     return;
4577   if ( e->_sWOL.IsNull() || e->_sWOL.ShapeType() != TopAbs_EDGE )
4578     throw SALOME_Exception(LOCALIZED("Wrong _LayerEdge is added"));
4579   if ( _edges[0] && _edges[0]->_sWOL != e->_sWOL )
4580     throw SALOME_Exception(LOCALIZED("Wrong _LayerEdge is added"));
4581
4582   // store _LayerEdge
4583   const TopoDS_Edge& E = TopoDS::Edge( e->_sWOL );
4584   double f,l;
4585   BRep_Tool::Range( E, f,l );
4586   double u = helper.GetNodeU( E, e->_nodes[0], e->_nodes.back());
4587   _edges[ u < 0.5*(f+l) ? 0 : 1 ] = e;
4588
4589   // Update _nodes
4590
4591   const SMDS_MeshNode* tgtNode0 = _edges[0] ? _edges[0]->_nodes.back() : 0;
4592   const SMDS_MeshNode* tgtNode1 = _edges[1] ? _edges[1]->_nodes.back() : 0;
4593
4594   if ( _nodes.empty() )
4595   {
4596     SMESHDS_SubMesh * eSubMesh = helper.GetMeshDS()->MeshElements( E );
4597     if ( !eSubMesh || eSubMesh->NbNodes() < 1 )
4598       return;
4599     TopLoc_Location loc;
4600     Handle(Geom_Curve) C = BRep_Tool::Curve(E, loc, f,l);
4601     GeomAdaptor_Curve aCurve(C, f,l);
4602     const double totLen = GCPnts_AbscissaPoint::Length(aCurve, f, l);
4603
4604     int nbExpectNodes = eSubMesh->NbNodes();
4605     _initU  .reserve( nbExpectNodes );
4606     _normPar.reserve( nbExpectNodes );
4607     _nodes  .reserve( nbExpectNodes );
4608     SMDS_NodeIteratorPtr nIt = eSubMesh->GetNodes();
4609     while ( nIt->more() )
4610     {
4611       const SMDS_MeshNode* node = nIt->next();
4612       if ( node->NbInverseElements(SMDSAbs_Edge) == 0 ||
4613            node == tgtNode0 || node == tgtNode1 )
4614         continue; // refinement nodes
4615       _nodes.push_back( node );
4616       _initU.push_back( helper.GetNodeU( E, node ));
4617       double len = GCPnts_AbscissaPoint::Length(aCurve, f, _initU.back());
4618       _normPar.push_back(  len / totLen );
4619     }
4620   }
4621   else
4622   {
4623     // remove target node of the _LayerEdge from _nodes
4624     int nbFound = 0;
4625     for ( size_t i = 0; i < _nodes.size(); ++i )
4626       if ( !_nodes[i] || _nodes[i] == tgtNode0 || _nodes[i] == tgtNode1 )
4627         _nodes[i] = 0, nbFound++;
4628     if ( nbFound == _nodes.size() )
4629       _nodes.clear();
4630   }
4631 }
4632
4633 //================================================================================
4634 /*!
4635  * \brief Move nodes on EDGE from ends where _LayerEdge's are inflated
4636  */
4637 //================================================================================
4638
4639 void _Shrinker1D::Compute(bool set3D, SMESH_MesherHelper& helper)
4640 {
4641   if ( _done || _nodes.empty())
4642     return;
4643   const _LayerEdge* e = _edges[0];
4644   if ( !e ) e = _edges[1];
4645   if ( !e ) return;
4646
4647   _done =  (( !_edges[0] || _edges[0]->_pos.empty() ) &&
4648             ( !_edges[1] || _edges[1]->_pos.empty() ));
4649
4650   const TopoDS_Edge& E = TopoDS::Edge( e->_sWOL );
4651   double f,l;
4652   if ( set3D || _done )
4653   {
4654     Handle(Geom_Curve) C = BRep_Tool::Curve(E, f,l);
4655     GeomAdaptor_Curve aCurve(C, f,l);
4656
4657     if ( _edges[0] )
4658       f = helper.GetNodeU( E, _edges[0]->_nodes.back(), _nodes[0] );
4659     if ( _edges[1] )
4660       l = helper.GetNodeU( E, _edges[1]->_nodes.back(), _nodes.back() );
4661     double totLen = GCPnts_AbscissaPoint::Length( aCurve, f, l );
4662
4663     for ( size_t i = 0; i < _nodes.size(); ++i )
4664     {
4665       if ( !_nodes[i] ) continue;
4666       double len = totLen * _normPar[i];
4667       GCPnts_AbscissaPoint discret( aCurve, len, f );
4668       if ( !discret.IsDone() )
4669         return throw SALOME_Exception(LOCALIZED("GCPnts_AbscissaPoint failed"));
4670       double u = discret.Parameter();
4671       SMDS_EdgePosition* pos = static_cast<SMDS_EdgePosition*>( _nodes[i]->GetPosition() );
4672       pos->SetUParameter( u );
4673       gp_Pnt p = C->Value( u );
4674       const_cast< SMDS_MeshNode*>( _nodes[i] )->setXYZ( p.X(), p.Y(), p.Z() );
4675     }
4676   }
4677   else
4678   {
4679     BRep_Tool::Range( E, f,l );
4680     if ( _edges[0] )
4681       f = helper.GetNodeU( E, _edges[0]->_nodes.back(), _nodes[0] );
4682     if ( _edges[1] )
4683       l = helper.GetNodeU( E, _edges[1]->_nodes.back(), _nodes.back() );
4684     
4685     for ( size_t i = 0; i < _nodes.size(); ++i )
4686     {
4687       if ( !_nodes[i] ) continue;
4688       double u = f * ( 1-_normPar[i] ) + l * _normPar[i];
4689       SMDS_EdgePosition* pos = static_cast<SMDS_EdgePosition*>( _nodes[i]->GetPosition() );
4690       pos->SetUParameter( u );
4691     }
4692   }
4693 }
4694
4695 //================================================================================
4696 /*!
4697  * \brief Restore initial parameters of nodes on EDGE
4698  */
4699 //================================================================================
4700
4701 void _Shrinker1D::RestoreParams()
4702 {
4703   if ( _done )
4704     for ( size_t i = 0; i < _nodes.size(); ++i )
4705     {
4706       if ( !_nodes[i] ) continue;
4707       SMDS_EdgePosition* pos = static_cast<SMDS_EdgePosition*>( _nodes[i]->GetPosition() );
4708       pos->SetUParameter( _initU[i] );
4709     }
4710   _done = false;
4711 }
4712
4713 //================================================================================
4714 /*!
4715  * \brief Replace source nodes by target nodes in shrinked mesh edges
4716  */
4717 //================================================================================
4718
4719 void _Shrinker1D::SwapSrcTgtNodes( SMESHDS_Mesh* mesh )
4720 {
4721   const SMDS_MeshNode* nodes[3];
4722   for ( int i = 0; i < 2; ++i )
4723   {
4724     if ( !_edges[i] ) continue;
4725
4726     SMESHDS_SubMesh * eSubMesh = mesh->MeshElements( _edges[i]->_sWOL );
4727     if ( !eSubMesh ) return;
4728     const SMDS_MeshNode* srcNode = _edges[i]->_nodes[0];
4729     const SMDS_MeshNode* tgtNode = _edges[i]->_nodes.back();
4730     SMDS_ElemIteratorPtr eIt = srcNode->GetInverseElementIterator(SMDSAbs_Edge);
4731     while ( eIt->more() )
4732     {
4733       const SMDS_MeshElement* e = eIt->next();
4734       if ( !eSubMesh->Contains( e ))
4735           continue;
4736       SMDS_ElemIteratorPtr nIt = e->nodesIterator();
4737       for ( int iN = 0; iN < e->NbNodes(); ++iN )
4738       {
4739         const SMDS_MeshNode* n = static_cast<const SMDS_MeshNode*>( nIt->next() );
4740         nodes[iN] = ( n == srcNode ? tgtNode : n );
4741       }
4742       mesh->ChangeElementNodes( e, nodes, e->NbNodes() );
4743     }
4744   }
4745 }
4746
4747 //================================================================================
4748 /*!
4749  * \brief Creates 2D and 1D elements on boundaries of new prisms
4750  */
4751 //================================================================================
4752
4753 bool _ViscousBuilder::addBoundaryElements()
4754 {
4755   SMESH_MesherHelper helper( *_mesh );
4756
4757   for ( size_t i = 0; i < _sdVec.size(); ++i )
4758   {
4759     _SolidData& data = _sdVec[i];
4760     TopTools_IndexedMapOfShape geomEdges;
4761     TopExp::MapShapes( data._solid, TopAbs_EDGE, geomEdges );
4762     for ( int iE = 1; iE <= geomEdges.Extent(); ++iE )
4763     {
4764       const TopoDS_Edge& E = TopoDS::Edge( geomEdges(iE));
4765
4766       // Get _LayerEdge's based on E
4767
4768       map< double, const SMDS_MeshNode* > u2nodes;
4769       if ( !SMESH_Algo::GetSortedNodesOnEdge( getMeshDS(), E, /*ignoreMedium=*/false, u2nodes))
4770         continue;
4771
4772       vector< _LayerEdge* > ledges; ledges.reserve( u2nodes.size() );
4773       TNode2Edge & n2eMap = data._n2eMap;
4774       map< double, const SMDS_MeshNode* >::iterator u2n = u2nodes.begin();
4775       {
4776         //check if 2D elements are needed on E
4777         TNode2Edge::iterator n2e = n2eMap.find( u2n->second );
4778         if ( n2e == n2eMap.end() ) continue; // no layers on vertex
4779         ledges.push_back( n2e->second );
4780         u2n++;
4781         if (( n2e = n2eMap.find( u2n->second )) == n2eMap.end() )
4782           continue; // no layers on E
4783         ledges.push_back( n2eMap[ u2n->second ]);
4784
4785         const SMDS_MeshNode* tgtN0 = ledges[0]->_nodes.back();
4786         const SMDS_MeshNode* tgtN1 = ledges[1]->_nodes.back();
4787         int nbSharedPyram = 0;
4788         SMDS_ElemIteratorPtr vIt = tgtN0->GetInverseElementIterator(SMDSAbs_Volume);
4789         while ( vIt->more() )
4790         {
4791           const SMDS_MeshElement* v = vIt->next();
4792           nbSharedPyram += int( v->GetNodeIndex( tgtN1 ) >= 0 );
4793         }
4794         if ( nbSharedPyram > 1 )
4795           continue; // not free border of the pyramid
4796
4797         if ( getMeshDS()->FindFace( ledges[0]->_nodes[0], ledges[0]->_nodes[1],
4798                                     ledges[1]->_nodes[0], ledges[1]->_nodes[1]))
4799           continue; // faces already created
4800       }
4801       for ( ++u2n; u2n != u2nodes.end(); ++u2n )
4802         ledges.push_back( n2eMap[ u2n->second ]);
4803
4804       // Find out orientation and type of face to create
4805
4806       bool reverse = false, isOnFace;
4807       
4808       map< TGeomID, TopoDS_Shape >::iterator e2f =
4809         data._shrinkShape2Shape.find( getMeshDS()->ShapeToIndex( E ));
4810       TopoDS_Shape F;
4811       if (( isOnFace = ( e2f != data._shrinkShape2Shape.end() )))
4812       {
4813         F = e2f->second.Oriented( TopAbs_FORWARD );
4814         reverse = ( helper.GetSubShapeOri( F, E ) == TopAbs_REVERSED );
4815         if ( helper.GetSubShapeOri( data._solid, F ) == TopAbs_REVERSED )
4816           reverse = !reverse, F.Reverse();
4817         if ( helper.IsReversedSubMesh( TopoDS::Face(F) ))
4818           reverse = !reverse;
4819       }
4820       else
4821       {
4822         // find FACE with layers sharing E
4823         PShapeIteratorPtr fIt = helper.GetAncestors( E, *_mesh, TopAbs_FACE );
4824         while ( fIt->more() && F.IsNull() )
4825         {
4826           const TopoDS_Shape* pF = fIt->next();
4827           if ( helper.IsSubShape( *pF, data._solid) &&
4828                !data._ignoreFaceIds.count( e2f->first ))
4829             F = *pF;
4830         }
4831       }
4832       // Find the sub-mesh to add new faces
4833       SMESHDS_SubMesh* sm = 0;
4834       if ( isOnFace )
4835         sm = getMeshDS()->MeshElements( F );
4836       else
4837         sm = data._proxyMesh->getFaceSubM( TopoDS::Face(F), /*create=*/true );
4838       if ( !sm )
4839         return error("error in addBoundaryElements()", data._index);
4840
4841       // Make faces
4842       const int dj1 = reverse ? 0 : 1;
4843       const int dj2 = reverse ? 1 : 0;
4844       for ( size_t j = 1; j < ledges.size(); ++j )
4845       {
4846         vector< const SMDS_MeshNode*>&  nn1 = ledges[j-dj1]->_nodes;
4847         vector< const SMDS_MeshNode*>&  nn2 = ledges[j-dj2]->_nodes;
4848         if ( isOnFace )
4849           for ( size_t z = 1; z < nn1.size(); ++z )
4850             sm->AddElement( getMeshDS()->AddFace( nn1[z-1], nn2[z-1], nn2[z], nn1[z] ));
4851         else
4852           for ( size_t z = 1; z < nn1.size(); ++z )
4853             sm->AddElement( new SMDS_FaceOfNodes( nn1[z-1], nn2[z-1], nn2[z], nn1[z]));
4854       }
4855
4856       // Make edges
4857       for ( int isFirst = 0; isFirst < 2; ++isFirst )
4858       {
4859         _LayerEdge* edge = isFirst ? ledges.front() : ledges.back();
4860         if ( !edge->_sWOL.IsNull() && edge->_sWOL.ShapeType() == TopAbs_EDGE )
4861         {
4862           vector< const SMDS_MeshNode*>&  nn = edge->_nodes;
4863           if ( nn[1]->GetInverseElementIterator( SMDSAbs_Edge )->more() )
4864             continue;
4865           helper.SetSubShape( edge->_sWOL );
4866           helper.SetElementsOnShape( true );
4867           for ( size_t z = 1; z < nn.size(); ++z )
4868             helper.AddEdge( nn[z-1], nn[z] );
4869         }
4870       }
4871     }
4872   }
4873
4874   return true;
4875 }