Salome HOME
#16522 [CEA 7599] Viscous layers hypothesis: extract layers as a group
[modules/smesh.git] / src / StdMeshers / StdMeshers_ViscousLayers.cxx
1 // Copyright (C) 2007-2019  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_PolygonalFaceOfNodes.hxx"
31 #include "SMDS_SetIterator.hxx"
32 #include "SMESHDS_Group.hxx"
33 #include "SMESHDS_Hypothesis.hxx"
34 #include "SMESHDS_Mesh.hxx"
35 #include "SMESH_Algo.hxx"
36 #include "SMESH_ComputeError.hxx"
37 #include "SMESH_ControlsDef.hxx"
38 #include "SMESH_Gen.hxx"
39 #include "SMESH_Group.hxx"
40 #include "SMESH_HypoFilter.hxx"
41 #include "SMESH_Mesh.hxx"
42 #include "SMESH_MeshAlgos.hxx"
43 #include "SMESH_MesherHelper.hxx"
44 #include "SMESH_ProxyMesh.hxx"
45 #include "SMESH_subMesh.hxx"
46 #include "SMESH_MeshEditor.hxx"
47 #include "SMESH_subMeshEventListener.hxx"
48 #include "StdMeshers_FaceSide.hxx"
49 #include "StdMeshers_ViscousLayers2D.hxx"
50
51 #include <Adaptor3d_HSurface.hxx>
52 #include <BRepAdaptor_Curve.hxx>
53 #include <BRepAdaptor_Curve2d.hxx>
54 #include <BRepAdaptor_Surface.hxx>
55 //#include <BRepLProp_CLProps.hxx>
56 #include <BRepLProp_SLProps.hxx>
57 #include <BRepOffsetAPI_MakeOffsetShape.hxx>
58 #include <BRep_Tool.hxx>
59 #include <Bnd_B2d.hxx>
60 #include <Bnd_B3d.hxx>
61 #include <ElCLib.hxx>
62 #include <GCPnts_AbscissaPoint.hxx>
63 #include <GCPnts_TangentialDeflection.hxx>
64 #include <Geom2d_Circle.hxx>
65 #include <Geom2d_Line.hxx>
66 #include <Geom2d_TrimmedCurve.hxx>
67 #include <GeomAdaptor_Curve.hxx>
68 #include <GeomLib.hxx>
69 #include <Geom_Circle.hxx>
70 #include <Geom_Curve.hxx>
71 #include <Geom_Line.hxx>
72 #include <Geom_TrimmedCurve.hxx>
73 #include <Precision.hxx>
74 #include <Standard_ErrorHandler.hxx>
75 #include <Standard_Failure.hxx>
76 #include <TColStd_Array1OfReal.hxx>
77 #include <TopExp.hxx>
78 #include <TopExp_Explorer.hxx>
79 #include <TopTools_IndexedMapOfShape.hxx>
80 #include <TopTools_ListOfShape.hxx>
81 #include <TopTools_MapIteratorOfMapOfShape.hxx>
82 #include <TopTools_MapOfShape.hxx>
83 #include <TopoDS.hxx>
84 #include <TopoDS_Edge.hxx>
85 #include <TopoDS_Face.hxx>
86 #include <TopoDS_Vertex.hxx>
87 #include <gp_Ax1.hxx>
88 #include <gp_Cone.hxx>
89 #include <gp_Sphere.hxx>
90 #include <gp_Vec.hxx>
91 #include <gp_XY.hxx>
92
93 #include <cmath>
94 #include <limits>
95 #include <list>
96 #include <queue>
97 #include <string>
98 #include <unordered_map>
99
100 #ifdef _DEBUG_
101 //#define __myDEBUG
102 //#define __NOT_INVALIDATE_BAD_SMOOTH
103 //#define __NODES_AT_POS
104 #endif
105
106 #define INCREMENTAL_SMOOTH // smooth only if min angle is too small
107 #define BLOCK_INFLATION // of individual _LayerEdge's
108 #define OLD_NEF_POLYGON
109
110 using namespace std;
111
112 //================================================================================
113 namespace VISCOUS_3D
114 {
115   typedef int TGeomID;
116
117   enum UIndex { U_TGT = 1, U_SRC, LEN_TGT };
118
119   const double theMinSmoothCosin = 0.1;
120   const double theSmoothThickToElemSizeRatio = 0.6;
121   const double theMinSmoothTriaAngle = 30;
122   const double theMinSmoothQuadAngle = 45;
123
124   // what part of thickness is allowed till intersection
125   // (defined by SALOME_TESTS/Grids/smesh/viscous_layers_00/A5)
126   const double theThickToIntersection = 1.5;
127
128   bool needSmoothing( double cosin, double tgtThick, double elemSize )
129   {
130     return cosin * tgtThick > theSmoothThickToElemSizeRatio * elemSize;
131   }
132   double getSmoothingThickness( double cosin, double elemSize )
133   {
134     return theSmoothThickToElemSizeRatio * elemSize / cosin;
135   }
136
137   /*!
138    * \brief SMESH_ProxyMesh computed by _ViscousBuilder for a SOLID.
139    * It is stored in a SMESH_subMesh of the SOLID as SMESH_subMeshEventListenerData
140    */
141   struct _MeshOfSolid : public SMESH_ProxyMesh,
142                         public SMESH_subMeshEventListenerData
143   {
144     bool                  _n2nMapComputed;
145     SMESH_ComputeErrorPtr _warning;
146
147     _MeshOfSolid( SMESH_Mesh* mesh)
148       :SMESH_subMeshEventListenerData( /*isDeletable=*/true),_n2nMapComputed(false)
149     {
150       SMESH_ProxyMesh::setMesh( *mesh );
151     }
152
153     // returns submesh for a geom face
154     SMESH_ProxyMesh::SubMesh* getFaceSubM(const TopoDS_Face& F, bool create=false)
155     {
156       TGeomID i = SMESH_ProxyMesh::shapeIndex(F);
157       return create ? SMESH_ProxyMesh::getProxySubMesh(i) : findProxySubMesh(i);
158     }
159     void setNode2Node(const SMDS_MeshNode*                 srcNode,
160                       const SMDS_MeshNode*                 proxyNode,
161                       const SMESH_ProxyMesh::SubMesh* subMesh)
162     {
163       SMESH_ProxyMesh::setNode2Node( srcNode,proxyNode,subMesh);
164     }
165   };
166   //--------------------------------------------------------------------------------
167   /*!
168    * \brief Listener of events of 3D sub-meshes computed with viscous layers.
169    * It is used to clear an inferior dim sub-meshes modified by viscous layers
170    */
171   class _ShrinkShapeListener : SMESH_subMeshEventListener
172   {
173     _ShrinkShapeListener()
174       : SMESH_subMeshEventListener(/*isDeletable=*/false,
175                                    "StdMeshers_ViscousLayers::_ShrinkShapeListener") {}
176   public:
177     static SMESH_subMeshEventListener* Get() { static _ShrinkShapeListener l; return &l; }
178     virtual void ProcessEvent(const int                       event,
179                               const int                       eventType,
180                               SMESH_subMesh*                  solidSM,
181                               SMESH_subMeshEventListenerData* data,
182                               const SMESH_Hypothesis*         hyp)
183     {
184       if ( SMESH_subMesh::COMPUTE_EVENT == eventType && solidSM->IsEmpty() && data )
185       {
186         SMESH_subMeshEventListener::ProcessEvent(event,eventType,solidSM,data,hyp);
187       }
188     }
189   };
190   //--------------------------------------------------------------------------------
191   /*!
192    * \brief Listener of events of 3D sub-meshes computed with viscous layers.
193    * It is used to store data computed by _ViscousBuilder for a sub-mesh and to
194    * delete the data as soon as it has been used
195    */
196   class _ViscousListener : SMESH_subMeshEventListener
197   {
198     _ViscousListener():
199       SMESH_subMeshEventListener(/*isDeletable=*/false,
200                                  "StdMeshers_ViscousLayers::_ViscousListener") {}
201     static SMESH_subMeshEventListener* Get() { static _ViscousListener l; return &l; }
202   public:
203     virtual void ProcessEvent(const int                       event,
204                               const int                       eventType,
205                               SMESH_subMesh*                  subMesh,
206                               SMESH_subMeshEventListenerData* data,
207                               const SMESH_Hypothesis*         hyp)
208     {
209       if (( SMESH_subMesh::COMPUTE_EVENT       == eventType ) &&
210           ( SMESH_subMesh::CHECK_COMPUTE_STATE != event &&
211             SMESH_subMesh::SUBMESH_COMPUTED    != event ))
212       {
213         // delete SMESH_ProxyMesh containing temporary faces
214         subMesh->DeleteEventListener( this );
215       }
216     }
217     // Finds or creates proxy mesh of the solid
218     static _MeshOfSolid* GetSolidMesh(SMESH_Mesh*         mesh,
219                                       const TopoDS_Shape& solid,
220                                       bool                toCreate=false)
221     {
222       if ( !mesh ) return 0;
223       SMESH_subMesh* sm = mesh->GetSubMesh(solid);
224       _MeshOfSolid* data = (_MeshOfSolid*) sm->GetEventListenerData( Get() );
225       if ( !data && toCreate )
226       {
227         data = new _MeshOfSolid(mesh);
228         data->mySubMeshes.push_back( sm ); // to find SOLID by _MeshOfSolid
229         sm->SetEventListener( Get(), data, sm );
230       }
231       return data;
232     }
233     // Removes proxy mesh of the solid
234     static void RemoveSolidMesh(SMESH_Mesh* mesh, const TopoDS_Shape& solid)
235     {
236       mesh->GetSubMesh(solid)->DeleteEventListener( _ViscousListener::Get() );
237     }
238   };
239   
240   //================================================================================
241   /*!
242    * \brief sets a sub-mesh event listener to clear sub-meshes of sub-shapes of
243    * the main shape when sub-mesh of the main shape is cleared,
244    * for example to clear sub-meshes of FACEs when sub-mesh of a SOLID
245    * is cleared
246    */
247   //================================================================================
248
249   void ToClearSubWithMain( SMESH_subMesh* sub, const TopoDS_Shape& main)
250   {
251     SMESH_subMesh* mainSM = sub->GetFather()->GetSubMesh( main );
252     SMESH_subMeshEventListenerData* data =
253       mainSM->GetEventListenerData( _ShrinkShapeListener::Get());
254     if ( data )
255     {
256       if ( find( data->mySubMeshes.begin(), data->mySubMeshes.end(), sub ) ==
257            data->mySubMeshes.end())
258         data->mySubMeshes.push_back( sub );
259     }
260     else
261     {
262       data = SMESH_subMeshEventListenerData::MakeData( /*dependent=*/sub );
263       sub->SetEventListener( _ShrinkShapeListener::Get(), data, /*whereToListenTo=*/mainSM );
264     }
265   }
266   struct _SolidData;
267   //--------------------------------------------------------------------------------
268   /*!
269    * \brief Simplex (triangle or tetrahedron) based on 1 (tria) or 2 (tet) nodes of
270    * _LayerEdge and 2 nodes of the mesh surface beening smoothed.
271    * The class is used to check validity of face or volumes around a smoothed node;
272    * it stores only 2 nodes as the other nodes are stored by _LayerEdge.
273    */
274   struct _Simplex
275   {
276     const SMDS_MeshNode *_nPrev, *_nNext; // nodes on a smoothed mesh surface
277     const SMDS_MeshNode *_nOpp; // in 2D case, a node opposite to a smoothed node in QUAD
278     _Simplex(const SMDS_MeshNode* nPrev=0,
279              const SMDS_MeshNode* nNext=0,
280              const SMDS_MeshNode* nOpp=0)
281       : _nPrev(nPrev), _nNext(nNext), _nOpp(nOpp) {}
282     bool IsForward(const gp_XYZ* pntSrc, const gp_XYZ* pntTgt, double& vol) const
283     {
284       const double M[3][3] =
285         {{ _nNext->X() - pntSrc->X(), _nNext->Y() - pntSrc->Y(), _nNext->Z() - pntSrc->Z() },
286          { pntTgt->X() - pntSrc->X(), pntTgt->Y() - pntSrc->Y(), pntTgt->Z() - pntSrc->Z() },
287          { _nPrev->X() - pntSrc->X(), _nPrev->Y() - pntSrc->Y(), _nPrev->Z() - pntSrc->Z() }};
288       vol = ( + M[0][0] * M[1][1] * M[2][2]
289               + M[0][1] * M[1][2] * M[2][0]
290               + M[0][2] * M[1][0] * M[2][1]
291               - M[0][0] * M[1][2] * M[2][1]
292               - M[0][1] * M[1][0] * M[2][2]
293               - M[0][2] * M[1][1] * M[2][0]);
294       return vol > 1e-100;
295     }
296     bool IsForward(const SMDS_MeshNode* nSrc, const gp_XYZ& pTgt, double& vol) const
297     {
298       SMESH_TNodeXYZ pSrc( nSrc );
299       return IsForward( &pSrc, &pTgt, vol );
300     }
301     bool IsForward(const gp_XY&         tgtUV,
302                    const SMDS_MeshNode* smoothedNode,
303                    const TopoDS_Face&   face,
304                    SMESH_MesherHelper&  helper,
305                    const double         refSign) const
306     {
307       gp_XY prevUV = helper.GetNodeUV( face, _nPrev, smoothedNode );
308       gp_XY nextUV = helper.GetNodeUV( face, _nNext, smoothedNode );
309       gp_Vec2d v1( tgtUV, prevUV ), v2( tgtUV, nextUV );
310       double d = v1 ^ v2;
311       return d*refSign > 1e-100;
312     }
313     bool IsMinAngleOK( const gp_XYZ& pTgt, double& minAngle ) const
314     {
315       SMESH_TNodeXYZ pPrev( _nPrev ), pNext( _nNext );
316       if ( !_nOpp ) // triangle
317       {
318         gp_Vec tp( pPrev - pTgt ), pn( pNext - pPrev ), nt( pTgt - pNext );
319         double tp2 = tp.SquareMagnitude();
320         double pn2 = pn.SquareMagnitude();
321         double nt2 = nt.SquareMagnitude();
322
323         if ( tp2 < pn2 && tp2 < nt2 )
324           minAngle = ( nt * -pn ) * ( nt * -pn ) / nt2 / pn2;
325         else if ( pn2 < nt2 )
326           minAngle = ( tp * -nt ) * ( tp * -nt ) / tp2 / nt2;
327         else
328           minAngle = ( pn * -tp ) * ( pn * -tp ) / pn2 / tp2;
329
330         static double theMaxCos2 = ( Cos( theMinSmoothTriaAngle * M_PI / 180. ) *
331                                      Cos( theMinSmoothTriaAngle * M_PI / 180. ));
332         return minAngle < theMaxCos2;
333       }
334       else // quadrangle
335       {
336         SMESH_TNodeXYZ pOpp( _nOpp );
337         gp_Vec tp( pPrev - pTgt ), po( pOpp - pPrev ), on( pNext - pOpp), nt( pTgt - pNext );
338         double tp2 = tp.SquareMagnitude();
339         double po2 = po.SquareMagnitude();
340         double on2 = on.SquareMagnitude();
341         double nt2 = nt.SquareMagnitude();
342         minAngle = Max( Max((( tp * -nt ) * ( tp * -nt ) / tp2 / nt2 ),
343                             (( po * -tp ) * ( po * -tp ) / po2 / tp2 )),
344                         Max((( on * -po ) * ( on * -po ) / on2 / po2 ),
345                             (( nt * -on ) * ( nt * -on ) / nt2 / on2 )));
346
347         static double theMaxCos2 = ( Cos( theMinSmoothQuadAngle * M_PI / 180. ) *
348                                      Cos( theMinSmoothQuadAngle * M_PI / 180. ));
349         return minAngle < theMaxCos2;
350       }
351     }
352     bool IsNeighbour(const _Simplex& other) const
353     {
354       return _nPrev == other._nNext || _nNext == other._nPrev;
355     }
356     bool Includes( const SMDS_MeshNode* node ) const { return _nPrev == node || _nNext == node; }
357     static void GetSimplices( const SMDS_MeshNode* node,
358                               vector<_Simplex>&   simplices,
359                               const set<TGeomID>& ingnoreShapes,
360                               const _SolidData*   dataToCheckOri = 0,
361                               const bool          toSort = false);
362     static void SortSimplices(vector<_Simplex>& simplices);
363   };
364   //--------------------------------------------------------------------------------
365   /*!
366    * Structure used to take into account surface curvature while smoothing
367    */
368   struct _Curvature
369   {
370     double   _r; // radius
371     double   _k; // factor to correct node smoothed position
372     double   _h2lenRatio; // avgNormProj / (2*avgDist)
373     gp_Pnt2d _uv; // UV used in putOnOffsetSurface()
374   public:
375     static _Curvature* New( double avgNormProj, double avgDist )
376     {
377       _Curvature* c = 0;
378       if ( fabs( avgNormProj / avgDist ) > 1./200 )
379       {
380         c = new _Curvature;
381         c->_r = avgDist * avgDist / avgNormProj;
382         c->_k = avgDist * avgDist / c->_r / c->_r;
383         //c->_k = avgNormProj / c->_r;
384         c->_k *= ( c->_r < 0 ? 1/1.1 : 1.1 ); // not to be too restrictive
385         c->_h2lenRatio = avgNormProj / ( avgDist + avgDist );
386
387         c->_uv.SetCoord( 0., 0. );
388       }
389       return c;
390     }
391     double lenDelta(double len) const { return _k * ( _r + len ); }
392     double lenDeltaByDist(double dist) const { return dist * _h2lenRatio; }
393   };
394   //--------------------------------------------------------------------------------
395
396   struct _2NearEdges;
397   struct _LayerEdge;
398   struct _EdgesOnShape;
399   struct _Smoother1D;
400   typedef map< const SMDS_MeshNode*, _LayerEdge*, TIDCompare > TNode2Edge;
401
402   //--------------------------------------------------------------------------------
403   /*!
404    * \brief Edge normal to surface, connecting a node on solid surface (_nodes[0])
405    * and a node of the most internal layer (_nodes.back())
406    */
407   struct _LayerEdge
408   {
409     typedef gp_XYZ (_LayerEdge::*PSmooFun)();
410
411     vector< const SMDS_MeshNode*> _nodes;
412
413     gp_XYZ              _normal;    // to boundary of solid
414     vector<gp_XYZ>      _pos;       // points computed during inflation
415     double              _len;       // length achieved with the last inflation step
416     double              _maxLen;    // maximal possible length
417     double              _cosin;     // of angle (_normal ^ surface)
418     double              _minAngle;  // of _simplices
419     double              _lenFactor; // to compute _len taking _cosin into account
420     int                 _flags;
421
422     // simplices connected to the source node (_nodes[0]);
423     // used for smoothing and quality check of _LayerEdge's based on the FACE
424     vector<_Simplex>    _simplices;
425     vector<_LayerEdge*> _neibors; // all surrounding _LayerEdge's
426     PSmooFun            _smooFunction; // smoothing function
427     _Curvature*         _curvature;
428     // data for smoothing of _LayerEdge's based on the EDGE
429     _2NearEdges*        _2neibors;
430
431     enum EFlags { TO_SMOOTH       = 0x0000001,
432                   MOVED           = 0x0000002, // set by _neibors[i]->SetNewLength()
433                   SMOOTHED        = 0x0000004, // set by _LayerEdge::Smooth()
434                   DIFFICULT       = 0x0000008, // near concave VERTEX
435                   ON_CONCAVE_FACE = 0x0000010,
436                   BLOCKED         = 0x0000020, // not to inflate any more
437                   INTERSECTED     = 0x0000040, // close intersection with a face found
438                   NORMAL_UPDATED  = 0x0000080,
439                   UPD_NORMAL_CONV = 0x0000100, // to update normal on boundary of concave FACE
440                   MARKED          = 0x0000200, // local usage
441                   MULTI_NORMAL    = 0x0000400, // a normal is invisible by some of surrounding faces
442                   NEAR_BOUNDARY   = 0x0000800, // is near FACE boundary forcing smooth
443                   SMOOTHED_C1     = 0x0001000, // is on _eosC1
444                   DISTORTED       = 0x0002000, // was bad before smoothing
445                   RISKY_SWOL      = 0x0004000, // SWOL is parallel to a source FACE
446                   SHRUNK          = 0x0008000, // target node reached a tgt position while shrink()
447                   UNUSED_FLAG     = 0x0100000  // to add user flags after
448     };
449     bool Is   ( int flag ) const { return _flags & flag; }
450     void Set  ( int flag ) { _flags |= flag; }
451     void Unset( int flag ) { _flags &= ~flag; }
452     std::string DumpFlags() const; // debug
453
454     void SetNewLength( double len, _EdgesOnShape& eos, SMESH_MesherHelper& helper );
455     bool SetNewLength2d( Handle(Geom_Surface)& surface,
456                          const TopoDS_Face&    F,
457                          _EdgesOnShape&        eos,
458                          SMESH_MesherHelper&   helper );
459     void SetDataByNeighbors( const SMDS_MeshNode* n1,
460                              const SMDS_MeshNode* n2,
461                              const _EdgesOnShape& eos,
462                              SMESH_MesherHelper&  helper);
463     void Block( _SolidData& data );
464     void InvalidateStep( size_t curStep, const _EdgesOnShape& eos, bool restoreLength=false );
465     void ChooseSmooFunction(const set< TGeomID >& concaveVertices,
466                             const TNode2Edge&     n2eMap);
467     void SmoothPos( const vector< double >& segLen, const double tol );
468     int  GetSmoothedPos( const double tol );
469     int  Smooth(const int step, const bool isConcaveFace, bool findBest);
470     int  Smooth(const int step, bool findBest, vector< _LayerEdge* >& toSmooth );
471     int  CheckNeiborsOnBoundary(vector< _LayerEdge* >* badNeibors = 0, bool * needSmooth = 0 );
472     void SmoothWoCheck();
473     bool SmoothOnEdge(Handle(ShapeAnalysis_Surface)& surface,
474                       const TopoDS_Face&             F,
475                       SMESH_MesherHelper&            helper);
476     void MoveNearConcaVer( const _EdgesOnShape*    eov,
477                            const _EdgesOnShape*    eos,
478                            const int               step,
479                            vector< _LayerEdge* > & badSmooEdges);
480     bool FindIntersection( SMESH_ElementSearcher&   searcher,
481                            double &                 distance,
482                            const double&            epsilon,
483                            _EdgesOnShape&           eos,
484                            const SMDS_MeshElement** face = 0);
485     bool SegTriaInter( const gp_Ax1&        lastSegment,
486                        const gp_XYZ&        p0,
487                        const gp_XYZ&        p1,
488                        const gp_XYZ&        p2,
489                        double&              dist,
490                        const double&        epsilon) const;
491     bool SegTriaInter( const gp_Ax1&        lastSegment,
492                        const SMDS_MeshNode* n0,
493                        const SMDS_MeshNode* n1,
494                        const SMDS_MeshNode* n2,
495                        double&              dist,
496                        const double&        epsilon) const
497     { return SegTriaInter( lastSegment,
498                            SMESH_TNodeXYZ( n0 ), SMESH_TNodeXYZ( n1 ), SMESH_TNodeXYZ( n2 ),
499                            dist, epsilon );
500     }
501     const gp_XYZ& PrevPos() const { return _pos[ _pos.size() - 2 ]; }
502     gp_XYZ PrevCheckPos( _EdgesOnShape* eos=0 ) const;
503     gp_Ax1 LastSegment(double& segLen, _EdgesOnShape& eos) const;
504     gp_XY  LastUV( const TopoDS_Face& F, _EdgesOnShape& eos, int which=-1 ) const;
505     bool   IsOnEdge() const { return _2neibors; }
506     bool   IsOnFace() const { return ( _nodes[0]->GetPosition()->GetDim() == 2 ); }
507     int    BaseShapeDim() const { return _nodes[0]->GetPosition()->GetDim(); }
508     gp_XYZ Copy( _LayerEdge& other, _EdgesOnShape& eos, SMESH_MesherHelper& helper );
509     void   SetCosin( double cosin );
510     void   SetNormal( const gp_XYZ& n ) { _normal = n; }
511     void   SetMaxLen( double l ) { _maxLen = l; }
512     int    NbSteps() const { return _pos.size() - 1; } // nb inlation steps
513     bool   IsNeiborOnEdge( const _LayerEdge* edge ) const;
514     void   SetSmooLen( double len ) { // set _len at which smoothing is needed
515       _cosin = len; // as for _LayerEdge's on FACE _cosin is not used
516     }
517     double GetSmooLen() { return _cosin; } // for _LayerEdge's on FACE _cosin is not used
518
519     gp_XYZ smoothLaplacian();
520     gp_XYZ smoothAngular();
521     gp_XYZ smoothLengthWeighted();
522     gp_XYZ smoothCentroidal();
523     gp_XYZ smoothNefPolygon();
524
525     enum { FUN_LAPLACIAN, FUN_LENWEIGHTED, FUN_CENTROIDAL, FUN_NEFPOLY, FUN_ANGULAR, FUN_NB };
526     static const int theNbSmooFuns = FUN_NB;
527     static PSmooFun _funs[theNbSmooFuns];
528     static const char* _funNames[theNbSmooFuns+1];
529     int smooFunID( PSmooFun fun=0) const;
530   };
531   _LayerEdge::PSmooFun _LayerEdge::_funs[theNbSmooFuns] = { &_LayerEdge::smoothLaplacian,
532                                                             &_LayerEdge::smoothLengthWeighted,
533                                                             &_LayerEdge::smoothCentroidal,
534                                                             &_LayerEdge::smoothNefPolygon,
535                                                             &_LayerEdge::smoothAngular };
536   const char* _LayerEdge::_funNames[theNbSmooFuns+1] = { "Laplacian",
537                                                          "LengthWeighted",
538                                                          "Centroidal",
539                                                          "NefPolygon",
540                                                          "Angular",
541                                                          "None"};
542   struct _LayerEdgeCmp
543   {
544     bool operator () (const _LayerEdge* e1, const _LayerEdge* e2) const
545     {
546       const bool cmpNodes = ( e1 && e2 && e1->_nodes.size() && e2->_nodes.size() );
547       return cmpNodes ? ( e1->_nodes[0]->GetID() < e2->_nodes[0]->GetID()) : ( e1 < e2 );
548     }
549   };
550   //--------------------------------------------------------------------------------
551   /*!
552    * A 2D half plane used by _LayerEdge::smoothNefPolygon()
553    */
554   struct _halfPlane
555   {
556     gp_XY _pos, _dir, _inNorm;
557     bool IsOut( const gp_XY p, const double tol ) const
558     {
559       return _inNorm * ( p - _pos ) < -tol;
560     }
561     bool FindIntersection( const _halfPlane& hp, gp_XY & intPnt )
562     {
563       //const double eps = 1e-10;
564       double D = _dir.Crossed( hp._dir );
565       if ( fabs(D) < std::numeric_limits<double>::min())
566         return false;
567       gp_XY vec21 = _pos - hp._pos; 
568       double u = hp._dir.Crossed( vec21 ) / D; 
569       intPnt = _pos + _dir * u;
570       return true;
571     }
572   };
573   //--------------------------------------------------------------------------------
574   /*!
575    * Structure used to smooth a _LayerEdge based on an EDGE.
576    */
577   struct _2NearEdges
578   {
579     double               _wgt  [2]; // weights of _nodes
580     _LayerEdge*          _edges[2];
581
582      // normal to plane passing through _LayerEdge._normal and tangent of EDGE
583     gp_XYZ*              _plnNorm;
584
585     _2NearEdges() { _edges[0]=_edges[1]=0; _plnNorm = 0; }
586     const SMDS_MeshNode* tgtNode(bool is2nd) {
587       return _edges[is2nd] ? _edges[is2nd]->_nodes.back() : 0;
588     }
589     const SMDS_MeshNode* srcNode(bool is2nd) {
590       return _edges[is2nd] ? _edges[is2nd]->_nodes[0] : 0;
591     }
592     void reverse() {
593       std::swap( _wgt  [0], _wgt  [1] );
594       std::swap( _edges[0], _edges[1] );
595     }
596     void set( _LayerEdge* e1, _LayerEdge* e2, double w1, double w2 ) {
597       _edges[0] = e1; _edges[1] = e2; _wgt[0] = w1; _wgt[1] = w2;
598     }
599     bool include( const _LayerEdge* e ) {
600       return ( _edges[0] == e || _edges[1] == e );
601     }
602   };
603
604
605   //--------------------------------------------------------------------------------
606   /*!
607    * \brief Layers parameters got by averaging several hypotheses
608    */
609   struct AverageHyp
610   {
611     AverageHyp( const StdMeshers_ViscousLayers* hyp = 0 )
612       :_nbLayers(0), _nbHyps(0), _method(0), _thickness(0), _stretchFactor(0)
613     {
614       Add( hyp );
615     }
616     void Add( const StdMeshers_ViscousLayers* hyp )
617     {
618       if ( hyp )
619       {
620         _nbHyps++;
621         _nbLayers       = hyp->GetNumberLayers();
622         //_thickness     += hyp->GetTotalThickness();
623         _thickness      = Max( _thickness, hyp->GetTotalThickness() );
624         _stretchFactor += hyp->GetStretchFactor();
625         _method         = hyp->GetMethod();
626         if ( _groupName.empty() )
627           _groupName = hyp->GetGroupName();
628       }
629     }
630     double GetTotalThickness() const { return _thickness; /*_nbHyps ? _thickness / _nbHyps : 0;*/ }
631     double GetStretchFactor()  const { return _nbHyps ? _stretchFactor / _nbHyps : 0; }
632     int    GetNumberLayers()   const { return _nbLayers; }
633     int    GetMethod()         const { return _method; }
634     bool   ToCreateGroup()     const { return !_groupName.empty(); }
635     const std::string& GetGroupName() const { return _groupName; }
636
637     bool   UseSurfaceNormal()  const
638     { return _method == StdMeshers_ViscousLayers::SURF_OFFSET_SMOOTH; }
639     bool   ToSmooth()          const
640     { return _method == StdMeshers_ViscousLayers::SURF_OFFSET_SMOOTH; }
641     bool   IsOffsetMethod()    const
642     { return _method == StdMeshers_ViscousLayers::FACE_OFFSET; }
643
644   private:
645     int         _nbLayers, _nbHyps, _method;
646     double      _thickness, _stretchFactor;
647     std::string _groupName;
648   };
649
650   //--------------------------------------------------------------------------------
651   /*!
652    * \brief _LayerEdge's on a shape and other shape data
653    */
654   struct _EdgesOnShape
655   {
656     vector< _LayerEdge* > _edges;
657
658     TopoDS_Shape          _shape;
659     TGeomID               _shapeID;
660     SMESH_subMesh *       _subMesh;
661     // face or edge w/o layer along or near which _edges are inflated
662     TopoDS_Shape          _sWOL;
663     bool                  _isRegularSWOL; // w/o singularities
664     // averaged StdMeshers_ViscousLayers parameters
665     AverageHyp            _hyp;
666     bool                  _toSmooth;
667     _Smoother1D*          _edgeSmoother;
668     vector< _EdgesOnShape* > _eosConcaVer; // edges at concave VERTEXes of a FACE
669     vector< _EdgesOnShape* > _eosC1; // to smooth together several C1 continues shapes
670
671     typedef std::unordered_map< const SMDS_MeshElement*, gp_XYZ > TFace2NormMap;
672     TFace2NormMap            _faceNormals; // if _shape is FACE
673     vector< _EdgesOnShape* > _faceEOS; // to get _faceNormals of adjacent FACEs
674
675     Handle(ShapeAnalysis_Surface) _offsetSurf;
676     _LayerEdge*                   _edgeForOffset;
677
678     _SolidData*            _data; // parent SOLID
679
680     _LayerEdge*      operator[](size_t i) const { return (_LayerEdge*) _edges[i]; }
681     size_t           size() const { return _edges.size(); }
682     TopAbs_ShapeEnum ShapeType() const
683     { return _shape.IsNull() ? TopAbs_SHAPE : _shape.ShapeType(); }
684     TopAbs_ShapeEnum SWOLType() const
685     { return _sWOL.IsNull() ? TopAbs_SHAPE : _sWOL.ShapeType(); }
686     bool             HasC1( const _EdgesOnShape* other ) const
687     { return std::find( _eosC1.begin(), _eosC1.end(), other ) != _eosC1.end(); }
688     bool             GetNormal( const SMDS_MeshElement* face, gp_Vec& norm );
689     _SolidData&      GetData() const { return *_data; }
690
691     _EdgesOnShape(): _shapeID(-1), _subMesh(0), _toSmooth(false), _edgeSmoother(0) {}
692   };
693
694   //--------------------------------------------------------------------------------
695   /*!
696    * \brief Convex FACE whose radius of curvature is less than the thickness of
697    *        layers. It is used to detect distortion of prisms based on a convex
698    *        FACE and to update normals to enable further increasing the thickness
699    */
700   struct _ConvexFace
701   {
702     TopoDS_Face                     _face;
703
704     // edges whose _simplices are used to detect prism distortion
705     vector< _LayerEdge* >           _simplexTestEdges;
706
707     // map a sub-shape to _SolidData::_edgesOnShape
708     map< TGeomID, _EdgesOnShape* >  _subIdToEOS;
709
710     bool                            _isTooCurved;
711     bool                            _normalsFixed;
712     bool                            _normalsFixedOnBorders; // used in putOnOffsetSurface()
713
714     double GetMaxCurvature( _SolidData&         data,
715                             _EdgesOnShape&      eof,
716                             BRepLProp_SLProps&  surfProp,
717                             SMESH_MesherHelper& helper);
718
719     bool GetCenterOfCurvature( _LayerEdge*         ledge,
720                                BRepLProp_SLProps&  surfProp,
721                                SMESH_MesherHelper& helper,
722                                gp_Pnt &            center ) const;
723     bool CheckPrisms() const;
724   };
725
726   //--------------------------------------------------------------------------------
727   /*!
728    * \brief Structure holding _LayerEdge's based on EDGEs that will collide
729    *        at inflation up to the full thickness. A detected collision
730    *        is fixed in updateNormals()
731    */
732   struct _CollisionEdges
733   {
734     _LayerEdge*           _edge;
735     vector< _LayerEdge* > _intEdges; // each pair forms an intersected quadrangle
736     const SMDS_MeshNode* nSrc(int i) const { return _intEdges[i]->_nodes[0]; }
737     const SMDS_MeshNode* nTgt(int i) const { return _intEdges[i]->_nodes.back(); }
738   };
739
740   //--------------------------------------------------------------------------------
741   /*!
742    * \brief Data of a SOLID
743    */
744   struct _SolidData
745   {
746     typedef const StdMeshers_ViscousLayers* THyp;
747     TopoDS_Shape                    _solid;
748     TopTools_MapOfShape             _before; // SOLIDs to be computed before _solid
749     TGeomID                         _index; // SOLID id
750     _MeshOfSolid*                   _proxyMesh;
751     list< THyp >                    _hyps;
752     list< TopoDS_Shape >            _hypShapes;
753     map< TGeomID, THyp >            _face2hyp; // filled if _hyps.size() > 1
754     set< TGeomID >                  _reversedFaceIds;
755     set< TGeomID >                  _ignoreFaceIds; // WOL FACEs and FACEs of other SOLIDs
756
757     double                          _stepSize, _stepSizeCoeff, _geomSize;
758     const SMDS_MeshNode*            _stepSizeNodes[2];
759
760     TNode2Edge                      _n2eMap; // nodes and _LayerEdge's based on them
761
762     // map to find _n2eMap of another _SolidData by a shrink shape shared by two _SolidData's
763     map< TGeomID, TNode2Edge* >     _s2neMap;
764     // _LayerEdge's with underlying shapes
765     vector< _EdgesOnShape >         _edgesOnShape;
766
767     // key:   an id of shape (EDGE or VERTEX) shared by a FACE with
768     //        layers and a FACE w/o layers
769     // value: the shape (FACE or EDGE) to shrink mesh on.
770     //       _LayerEdge's basing on nodes on key shape are inflated along the value shape
771     map< TGeomID, TopoDS_Shape >     _shrinkShape2Shape;
772
773     // Convex FACEs whose radius of curvature is less than the thickness of layers
774     map< TGeomID, _ConvexFace >      _convexFaces;
775
776     // shapes (EDGEs and VERTEXes) shrink from which is forbidden due to collisions with
777     // the adjacent SOLID
778     set< TGeomID >                   _noShrinkShapes;
779
780     int                              _nbShapesToSmooth;
781
782     vector< _CollisionEdges >        _collisionEdges;
783     set< TGeomID >                   _concaveFaces;
784
785     double                           _maxThickness; // of all _hyps
786     double                           _minThickness; // of all _hyps
787
788     double                           _epsilon; // precision for SegTriaInter()
789
790     SMESH_MesherHelper*              _helper;
791
792     _SolidData(const TopoDS_Shape& s=TopoDS_Shape(),
793                _MeshOfSolid*       m=0)
794       :_solid(s), _proxyMesh(m), _helper(0) {}
795     ~_SolidData();
796
797     void SortOnEdge( const TopoDS_Edge& E, vector< _LayerEdge* >& edges);
798     void Sort2NeiborsOnEdge( vector< _LayerEdge* >& edges );
799
800     _ConvexFace* GetConvexFace( const TGeomID faceID ) {
801       map< TGeomID, _ConvexFace >::iterator id2face = _convexFaces.find( faceID );
802       return id2face == _convexFaces.end() ? 0 : & id2face->second;
803     }
804     _EdgesOnShape* GetShapeEdges(const TGeomID       shapeID );
805     _EdgesOnShape* GetShapeEdges(const TopoDS_Shape& shape );
806     _EdgesOnShape* GetShapeEdges(const _LayerEdge*   edge )
807     { return GetShapeEdges( edge->_nodes[0]->getshapeId() ); }
808
809     SMESH_MesherHelper& GetHelper() const { return *_helper; }
810
811     void UnmarkEdges( int flag = _LayerEdge::MARKED ) {
812       for ( size_t i = 0; i < _edgesOnShape.size(); ++i )
813         for ( size_t j = 0; j < _edgesOnShape[i]._edges.size(); ++j )
814           _edgesOnShape[i]._edges[j]->Unset( flag );
815     }
816     void AddShapesToSmooth( const set< _EdgesOnShape* >& shape,
817                             const set< _EdgesOnShape* >* edgesNoAnaSmooth=0 );
818
819     void PrepareEdgesToSmoothOnFace( _EdgesOnShape* eof, bool substituteSrcNodes );
820   };
821   //--------------------------------------------------------------------------------
822   /*!
823    * \brief Offset plane used in getNormalByOffset()
824    */
825   struct _OffsetPlane
826   {
827     gp_Pln _plane;
828     int    _faceIndex;
829     int    _faceIndexNext[2];
830     gp_Lin _lines[2]; // line of intersection with neighbor _OffsetPlane's
831     bool   _isLineOK[2];
832     _OffsetPlane() {
833       _isLineOK[0] = _isLineOK[1] = false; _faceIndexNext[0] = _faceIndexNext[1] = -1;
834     }
835     void   ComputeIntersectionLine( _OffsetPlane&        pln, 
836                                     const TopoDS_Edge&   E,
837                                     const TopoDS_Vertex& V );
838     gp_XYZ GetCommonPoint(bool& isFound, const TopoDS_Vertex& V) const;
839     int    NbLines() const { return _isLineOK[0] + _isLineOK[1]; }
840   };
841   //--------------------------------------------------------------------------------
842   /*!
843    * \brief Container of centers of curvature at nodes on an EDGE bounding _ConvexFace
844    */
845   struct _CentralCurveOnEdge
846   {
847     bool                  _isDegenerated;
848     vector< gp_Pnt >      _curvaCenters;
849     vector< _LayerEdge* > _ledges;
850     vector< gp_XYZ >      _normals; // new normal for each of _ledges
851     vector< double >      _segLength2;
852
853     TopoDS_Edge           _edge;
854     TopoDS_Face           _adjFace;
855     bool                  _adjFaceToSmooth;
856
857     void Append( const gp_Pnt& center, _LayerEdge* ledge )
858     {
859       if ( ledge->Is( _LayerEdge::MULTI_NORMAL ))
860         return;
861       if ( _curvaCenters.size() > 0 )
862         _segLength2.push_back( center.SquareDistance( _curvaCenters.back() ));
863       _curvaCenters.push_back( center );
864       _ledges.push_back( ledge );
865       _normals.push_back( ledge->_normal );
866     }
867     bool FindNewNormal( const gp_Pnt& center, gp_XYZ& newNormal );
868     void SetShapes( const TopoDS_Edge&  edge,
869                     const _ConvexFace&  convFace,
870                     _SolidData&         data,
871                     SMESH_MesherHelper& helper);
872   };
873   //--------------------------------------------------------------------------------
874   /*!
875    * \brief Data of node on a shrinked FACE
876    */
877   struct _SmoothNode
878   {
879     const SMDS_MeshNode*         _node;
880     vector<_Simplex>             _simplices; // for quality check
881
882     enum SmoothType { LAPLACIAN, CENTROIDAL, ANGULAR, TFI };
883
884     bool Smooth(int&                  badNb,
885                 Handle(Geom_Surface)& surface,
886                 SMESH_MesherHelper&   helper,
887                 const double          refSign,
888                 SmoothType            how,
889                 bool                  set3D);
890
891     gp_XY computeAngularPos(vector<gp_XY>& uv,
892                             const gp_XY&   uvToFix,
893                             const double   refSign );
894   };
895   struct PyDump;
896   //--------------------------------------------------------------------------------
897   /*!
898    * \brief Builder of viscous layers
899    */
900   class _ViscousBuilder
901   {
902   public:
903     _ViscousBuilder();
904     // does it's job
905     SMESH_ComputeErrorPtr Compute(SMESH_Mesh&         mesh,
906                                   const TopoDS_Shape& shape);
907     // check validity of hypotheses
908     SMESH_ComputeErrorPtr CheckHypotheses( SMESH_Mesh&         mesh,
909                                            const TopoDS_Shape& shape );
910
911     // restore event listeners used to clear an inferior dim sub-mesh modified by viscous layers
912     void RestoreListeners();
913
914     // computes SMESH_ProxyMesh::SubMesh::_n2n;
915     bool MakeN2NMap( _MeshOfSolid* pm );
916
917   private:
918
919     bool findSolidsWithLayers();
920     bool setBefore( _SolidData& solidBefore, _SolidData& solidAfter );
921     bool findFacesWithLayers(const bool onlyWith=false);
922     void getIgnoreFaces(const TopoDS_Shape&             solid,
923                         const StdMeshers_ViscousLayers* hyp,
924                         const TopoDS_Shape&             hypShape,
925                         set<TGeomID>&                   ignoreFaces);
926     bool makeLayer(_SolidData& data);
927     void setShapeData( _EdgesOnShape& eos, SMESH_subMesh* sm, _SolidData& data );
928     bool setEdgeData( _LayerEdge& edge, _EdgesOnShape& eos,
929                       SMESH_MesherHelper& helper, _SolidData& data);
930     gp_XYZ getFaceNormal(const SMDS_MeshNode* n,
931                          const TopoDS_Face&   face,
932                          SMESH_MesherHelper&  helper,
933                          bool&                isOK,
934                          bool                 shiftInside=false);
935     bool getFaceNormalAtSingularity(const gp_XY&        uv,
936                                     const TopoDS_Face&  face,
937                                     SMESH_MesherHelper& helper,
938                                     gp_Dir&             normal );
939     gp_XYZ getWeigthedNormal( const _LayerEdge*                edge );
940     gp_XYZ getNormalByOffset( _LayerEdge*                      edge,
941                               std::pair< TopoDS_Face, gp_XYZ > fId2Normal[],
942                               int                              nbFaces,
943                               bool                             lastNoOffset = false);
944     bool findNeiborsOnEdge(const _LayerEdge*     edge,
945                            const SMDS_MeshNode*& n1,
946                            const SMDS_MeshNode*& n2,
947                            _EdgesOnShape&        eos,
948                            _SolidData&           data);
949     void findSimplexTestEdges( _SolidData&                    data,
950                                vector< vector<_LayerEdge*> >& edgesByGeom);
951     void computeGeomSize( _SolidData& data );
952     bool findShapesToSmooth( _SolidData& data);
953     void limitStepSizeByCurvature( _SolidData&  data );
954     void limitStepSize( _SolidData&             data,
955                         const SMDS_MeshElement* face,
956                         const _LayerEdge*       maxCosinEdge );
957     void limitStepSize( _SolidData& data, const double minSize);
958     bool inflate(_SolidData& data);
959     bool smoothAndCheck(_SolidData& data, const int nbSteps, double & distToIntersection);
960     int  invalidateBadSmooth( _SolidData&               data,
961                               SMESH_MesherHelper&       helper,
962                               vector< _LayerEdge* >&    badSmooEdges,
963                               vector< _EdgesOnShape* >& eosC1,
964                               const int                 infStep );
965     void makeOffsetSurface( _EdgesOnShape& eos, SMESH_MesherHelper& );
966     void putOnOffsetSurface( _EdgesOnShape& eos, int infStep,
967                              vector< _EdgesOnShape* >& eosC1,
968                              int smooStep=0, int moveAll=false );
969     void findCollisionEdges( _SolidData& data, SMESH_MesherHelper& helper );
970     void findEdgesToUpdateNormalNearConvexFace( _ConvexFace &       convFace,
971                                                 _SolidData&         data,
972                                                 SMESH_MesherHelper& helper );
973     void limitMaxLenByCurvature( _SolidData& data, SMESH_MesherHelper& helper );
974     void limitMaxLenByCurvature( _LayerEdge* e1, _LayerEdge* e2,
975                                  _EdgesOnShape& eos1, _EdgesOnShape& eos2,
976                                  const bool isSmoothable );
977     bool updateNormals( _SolidData& data, SMESH_MesherHelper& helper, int stepNb, double stepSize );
978     bool updateNormalsOfConvexFaces( _SolidData&         data,
979                                      SMESH_MesherHelper& helper,
980                                      int                 stepNb );
981     void updateNormalsOfC1Vertices( _SolidData& data );
982     bool updateNormalsOfSmoothed( _SolidData&         data,
983                                   SMESH_MesherHelper& helper,
984                                   const int           nbSteps,
985                                   const double        stepSize );
986     bool isNewNormalOk( _SolidData&   data,
987                         _LayerEdge&   edge,
988                         const gp_XYZ& newNormal);
989     bool refine(_SolidData& data);
990     bool shrink(_SolidData& data);
991     bool prepareEdgeToShrink( _LayerEdge& edge, _EdgesOnShape& eos,
992                               SMESH_MesherHelper& helper,
993                               const SMESHDS_SubMesh* faceSubMesh );
994     void restoreNoShrink( _LayerEdge& edge ) const;
995     void fixBadFaces(const TopoDS_Face&          F,
996                      SMESH_MesherHelper&         helper,
997                      const bool                  is2D,
998                      const int                   step,
999                      set<const SMDS_MeshNode*> * involvedNodes=NULL);
1000     bool addBoundaryElements(_SolidData& data);
1001
1002     bool error( const string& text, int solidID=-1 );
1003     SMESHDS_Mesh* getMeshDS() const { return _mesh->GetMeshDS(); }
1004
1005     // debug
1006     void makeGroupOfLE();
1007
1008     SMESH_Mesh*                _mesh;
1009     SMESH_ComputeErrorPtr      _error;
1010
1011     vector<                    _SolidData >  _sdVec;
1012     TopTools_IndexedMapOfShape _solids; // to find _SolidData by a solid
1013     TopTools_MapOfShape        _shrinkedFaces;
1014
1015     int                        _tmpFaceID;
1016     PyDump*                    _pyDump;
1017   };
1018   //--------------------------------------------------------------------------------
1019   /*!
1020    * \brief Shrinker of nodes on the EDGE
1021    */
1022   class _Shrinker1D
1023   {
1024     TopoDS_Edge                   _geomEdge;
1025     vector<double>                _initU;
1026     vector<double>                _normPar;
1027     vector<const SMDS_MeshNode*>  _nodes;
1028     const _LayerEdge*             _edges[2];
1029     bool                          _done;
1030   public:
1031     void AddEdge( const _LayerEdge* e, _EdgesOnShape& eos, SMESH_MesherHelper& helper );
1032     void Compute(bool set3D, SMESH_MesherHelper& helper);
1033     void RestoreParams();
1034     void SwapSrcTgtNodes(SMESHDS_Mesh* mesh);
1035     const TopoDS_Edge& GeomEdge() const { return _geomEdge; }
1036     const SMDS_MeshNode* TgtNode( bool is2nd ) const
1037     { return _edges[is2nd] ? _edges[is2nd]->_nodes.back() : 0; }
1038     const SMDS_MeshNode* SrcNode( bool is2nd ) const
1039     { return _edges[is2nd] ? _edges[is2nd]->_nodes[0] : 0; }
1040   };
1041   //--------------------------------------------------------------------------------
1042   /*!
1043    * \brief Smoother of _LayerEdge's on EDGE.
1044    */
1045   struct _Smoother1D
1046   {
1047     struct OffPnt // point of the offsetted EDGE
1048     {
1049       gp_XYZ      _xyz;    // coord of a point inflated from EDGE w/o smooth
1050       double      _len;    // length reached at previous inflation step
1051       double      _param;  // on EDGE
1052       _2NearEdges _2edges; // 2 neighbor _LayerEdge's
1053       gp_XYZ      _edgeDir;// EDGE tangent at _param
1054       double Distance( const OffPnt& p ) const { return ( _xyz - p._xyz ).Modulus(); }
1055     };
1056     vector< OffPnt >   _offPoints;
1057     vector< double >   _leParams; // normalized param of _eos._edges on EDGE
1058     Handle(Geom_Curve) _anaCurve; // for analytic smooth
1059     _LayerEdge         _leOnV[2]; // _LayerEdge's holding normal to the EDGE at VERTEXes
1060     gp_XYZ             _edgeDir[2]; // tangent at VERTEXes
1061     size_t             _iSeg[2];  // index of segment where extreme tgt node is projected
1062     _EdgesOnShape&     _eos;
1063     double             _curveLen; // length of the EDGE
1064     std::pair<int,int> _eToSmooth[2]; // <from,to> indices of _LayerEdge's in _eos
1065
1066     static Handle(Geom_Curve) CurveForSmooth( const TopoDS_Edge&  E,
1067                                               _EdgesOnShape&      eos,
1068                                               SMESH_MesherHelper& helper);
1069
1070     _Smoother1D( Handle(Geom_Curve) curveForSmooth,
1071                  _EdgesOnShape&     eos )
1072       : _anaCurve( curveForSmooth ), _eos( eos )
1073     {
1074     }
1075     bool Perform(_SolidData&                    data,
1076                  Handle(ShapeAnalysis_Surface)& surface,
1077                  const TopoDS_Face&             F,
1078                  SMESH_MesherHelper&            helper );
1079
1080     void prepare(_SolidData& data );
1081
1082     void findEdgesToSmooth();
1083
1084     bool isToSmooth( int iE );
1085
1086     bool smoothAnalyticEdge( _SolidData&                    data,
1087                              Handle(ShapeAnalysis_Surface)& surface,
1088                              const TopoDS_Face&             F,
1089                              SMESH_MesherHelper&            helper);
1090     bool smoothComplexEdge( _SolidData&                    data,
1091                             Handle(ShapeAnalysis_Surface)& surface,
1092                             const TopoDS_Face&             F,
1093                             SMESH_MesherHelper&            helper);
1094     gp_XYZ getNormalNormal( const gp_XYZ & normal,
1095                             const gp_XYZ&  edgeDir);
1096     _LayerEdge* getLEdgeOnV( bool is2nd )
1097     {
1098       return _eos._edges[ is2nd ? _eos._edges.size()-1 : 0 ]->_2neibors->_edges[ is2nd ];
1099     }
1100     bool isAnalytic() const { return !_anaCurve.IsNull(); }
1101
1102     void offPointsToPython() const; // debug
1103   };
1104   //--------------------------------------------------------------------------------
1105   /*!
1106    * \brief Class of temporary mesh face.
1107    * We can't use SMDS_FaceOfNodes since it's impossible to set it's ID which is
1108    * needed because SMESH_ElementSearcher internally uses set of elements sorted by ID
1109    */
1110   struct _TmpMeshFace : public SMDS_PolygonalFaceOfNodes
1111   {
1112     const SMDS_MeshElement* _srcFace;
1113
1114     _TmpMeshFace( const vector<const SMDS_MeshNode*>& nodes,
1115                   int                                 ID,
1116                   int                                 faceID=-1,
1117                   const SMDS_MeshElement*             srcFace=0 ):
1118       SMDS_PolygonalFaceOfNodes(nodes), _srcFace( srcFace ) { setID( ID ); setShapeID( faceID ); }
1119     virtual SMDSAbs_EntityType  GetEntityType() const
1120     { return _srcFace ? _srcFace->GetEntityType() : SMDSEntity_Quadrangle; }
1121     virtual SMDSAbs_GeometryType GetGeomType()  const
1122     { return _srcFace ? _srcFace->GetGeomType() : SMDSGeom_QUADRANGLE; }
1123   };
1124   //--------------------------------------------------------------------------------
1125   /*!
1126    * \brief Class of temporary mesh quadrangle face storing _LayerEdge it's based on
1127    */
1128   struct _TmpMeshFaceOnEdge : public _TmpMeshFace
1129   {
1130     _LayerEdge *_le1, *_le2;
1131     _TmpMeshFaceOnEdge( _LayerEdge* le1, _LayerEdge* le2, int ID ):
1132       _TmpMeshFace( vector<const SMDS_MeshNode*>(4), ID ), _le1(le1), _le2(le2)
1133     {
1134       myNodes[0]=_le1->_nodes[0];
1135       myNodes[1]=_le1->_nodes.back();
1136       myNodes[2]=_le2->_nodes.back();
1137       myNodes[3]=_le2->_nodes[0];
1138     }
1139     const SMDS_MeshNode* n( size_t i ) const
1140     {
1141       return myNodes[ i ];
1142     }
1143     gp_XYZ GetDir() const // return average direction of _LayerEdge's, normal to EDGE
1144     {
1145       SMESH_TNodeXYZ p0s( myNodes[0] );
1146       SMESH_TNodeXYZ p0t( myNodes[1] );
1147       SMESH_TNodeXYZ p1t( myNodes[2] );
1148       SMESH_TNodeXYZ p1s( myNodes[3] );
1149       gp_XYZ  v0 = p0t - p0s;
1150       gp_XYZ  v1 = p1t - p1s;
1151       gp_XYZ v01 = p1s - p0s;
1152       gp_XYZ   n = ( v0 ^ v01 ) + ( v1 ^ v01 );
1153       gp_XYZ   d = v01 ^ n;
1154       d.Normalize();
1155       return d;
1156     }
1157     gp_XYZ GetDir(_LayerEdge* le1, _LayerEdge* le2) // return average direction of _LayerEdge's
1158     {
1159       myNodes[0]=le1->_nodes[0];
1160       myNodes[1]=le1->_nodes.back();
1161       myNodes[2]=le2->_nodes.back();
1162       myNodes[3]=le2->_nodes[0];
1163       return GetDir();
1164     }
1165   };
1166   //--------------------------------------------------------------------------------
1167   /*!
1168    * \brief Retriever of node coordinates either directly or from a surface by node UV.
1169    * \warning Location of a surface is ignored
1170    */
1171   struct _NodeCoordHelper
1172   {
1173     SMESH_MesherHelper&        _helper;
1174     const TopoDS_Face&         _face;
1175     Handle(Geom_Surface)       _surface;
1176     gp_XYZ (_NodeCoordHelper::* _fun)(const SMDS_MeshNode* n) const;
1177
1178     _NodeCoordHelper(const TopoDS_Face& F, SMESH_MesherHelper& helper, bool is2D)
1179       : _helper( helper ), _face( F )
1180     {
1181       if ( is2D )
1182       {
1183         TopLoc_Location loc;
1184         _surface = BRep_Tool::Surface( _face, loc );
1185       }
1186       if ( _surface.IsNull() )
1187         _fun = & _NodeCoordHelper::direct;
1188       else
1189         _fun = & _NodeCoordHelper::byUV;
1190     }
1191     gp_XYZ operator()(const SMDS_MeshNode* n) const { return (this->*_fun)( n ); }
1192
1193   private:
1194     gp_XYZ direct(const SMDS_MeshNode* n) const
1195     {
1196       return SMESH_TNodeXYZ( n );
1197     }
1198     gp_XYZ byUV  (const SMDS_MeshNode* n) const
1199     {
1200       gp_XY uv = _helper.GetNodeUV( _face, n );
1201       return _surface->Value( uv.X(), uv.Y() ).XYZ();
1202     }
1203   };
1204
1205   //================================================================================
1206   /*!
1207    * \brief Check angle between vectors 
1208    */
1209   //================================================================================
1210
1211   inline bool isLessAngle( const gp_Vec& v1, const gp_Vec& v2, const double cos )
1212   {
1213     double dot = v1 * v2; // cos * |v1| * |v2|
1214     double l1  = v1.SquareMagnitude();
1215     double l2  = v2.SquareMagnitude();
1216     return (( dot * cos >= 0 ) && 
1217             ( dot * dot ) / l1 / l2 >= ( cos * cos ));
1218   }
1219
1220 } // namespace VISCOUS_3D
1221
1222
1223
1224 //================================================================================
1225 // StdMeshers_ViscousLayers hypothesis
1226 //
1227 StdMeshers_ViscousLayers::StdMeshers_ViscousLayers(int hypId, SMESH_Gen* gen)
1228   :SMESH_Hypothesis(hypId, gen),
1229    _isToIgnoreShapes(1), _nbLayers(1), _thickness(1), _stretchFactor(1),
1230    _method( SURF_OFFSET_SMOOTH ),
1231    _groupName("")
1232 {
1233   _name = StdMeshers_ViscousLayers::GetHypType();
1234   _param_algo_dim = -3; // auxiliary hyp used by 3D algos
1235 } // --------------------------------------------------------------------------------
1236 void StdMeshers_ViscousLayers::SetBndShapes(const std::vector<int>& faceIds, bool toIgnore)
1237 {
1238   if ( faceIds != _shapeIds )
1239     _shapeIds = faceIds, NotifySubMeshesHypothesisModification();
1240   if ( _isToIgnoreShapes != toIgnore )
1241     _isToIgnoreShapes = toIgnore, NotifySubMeshesHypothesisModification();
1242 } // --------------------------------------------------------------------------------
1243 void StdMeshers_ViscousLayers::SetTotalThickness(double thickness)
1244 {
1245   if ( thickness != _thickness )
1246     _thickness = thickness, NotifySubMeshesHypothesisModification();
1247 } // --------------------------------------------------------------------------------
1248 void StdMeshers_ViscousLayers::SetNumberLayers(int nb)
1249 {
1250   if ( _nbLayers != nb )
1251     _nbLayers = nb, NotifySubMeshesHypothesisModification();
1252 } // --------------------------------------------------------------------------------
1253 void StdMeshers_ViscousLayers::SetStretchFactor(double factor)
1254 {
1255   if ( _stretchFactor != factor )
1256     _stretchFactor = factor, NotifySubMeshesHypothesisModification();
1257 } // --------------------------------------------------------------------------------
1258 void StdMeshers_ViscousLayers::SetMethod( ExtrusionMethod method )
1259 {
1260   if ( _method != method )
1261     _method = method, NotifySubMeshesHypothesisModification();
1262 } // --------------------------------------------------------------------------------
1263 void StdMeshers_ViscousLayers::SetGroupName(const std::string& name)
1264 {
1265   if ( _groupName != name )
1266   {
1267     _groupName = name;
1268     if ( !_groupName.empty() )
1269       NotifySubMeshesHypothesisModification();
1270   }
1271 } // --------------------------------------------------------------------------------
1272 SMESH_ProxyMesh::Ptr
1273 StdMeshers_ViscousLayers::Compute(SMESH_Mesh&         theMesh,
1274                                   const TopoDS_Shape& theShape,
1275                                   const bool          toMakeN2NMap) const
1276 {
1277   using namespace VISCOUS_3D;
1278   _ViscousBuilder builder;
1279   SMESH_ComputeErrorPtr err = builder.Compute( theMesh, theShape );
1280   if ( err && !err->IsOK() )
1281     return SMESH_ProxyMesh::Ptr();
1282
1283   vector<SMESH_ProxyMesh::Ptr> components;
1284   TopExp_Explorer exp( theShape, TopAbs_SOLID );
1285   for ( ; exp.More(); exp.Next() )
1286   {
1287     if ( _MeshOfSolid* pm =
1288          _ViscousListener::GetSolidMesh( &theMesh, exp.Current(), /*toCreate=*/false))
1289     {
1290       if ( toMakeN2NMap && !pm->_n2nMapComputed )
1291         if ( !builder.MakeN2NMap( pm ))
1292           return SMESH_ProxyMesh::Ptr();
1293       components.push_back( SMESH_ProxyMesh::Ptr( pm ));
1294       pm->myIsDeletable = false; // it will de deleted by boost::shared_ptr
1295
1296       if ( pm->_warning && !pm->_warning->IsOK() )
1297       {
1298         SMESH_subMesh* sm = theMesh.GetSubMesh( exp.Current() );
1299         SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1300         if ( !smError || smError->IsOK() )
1301           smError = pm->_warning;
1302       }
1303     }
1304     _ViscousListener::RemoveSolidMesh ( &theMesh, exp.Current() );
1305   }
1306   switch ( components.size() )
1307   {
1308   case 0: break;
1309
1310   case 1: return components[0];
1311
1312   default: return SMESH_ProxyMesh::Ptr( new SMESH_ProxyMesh( components ));
1313   }
1314   return SMESH_ProxyMesh::Ptr();
1315 } // --------------------------------------------------------------------------------
1316 std::ostream & StdMeshers_ViscousLayers::SaveTo(std::ostream & save)
1317 {
1318   save << " " << _nbLayers
1319        << " " << _thickness
1320        << " " << _stretchFactor
1321        << " " << _shapeIds.size();
1322   for ( size_t i = 0; i < _shapeIds.size(); ++i )
1323     save << " " << _shapeIds[i];
1324   save << " " << !_isToIgnoreShapes; // negate to keep the behavior in old studies.
1325   save << " " << _method;
1326   save << " " << _groupName.size();
1327   if ( !_groupName.empty() )
1328     save << " " << _groupName;
1329   return save;
1330 } // --------------------------------------------------------------------------------
1331 std::istream & StdMeshers_ViscousLayers::LoadFrom(std::istream & load)
1332 {
1333   int nbFaces, faceID, shapeToTreat, method;
1334   load >> _nbLayers >> _thickness >> _stretchFactor >> nbFaces;
1335   while ( (int) _shapeIds.size() < nbFaces && load >> faceID )
1336     _shapeIds.push_back( faceID );
1337   if ( load >> shapeToTreat ) {
1338     _isToIgnoreShapes = !shapeToTreat;
1339     if ( load >> method )
1340       _method = (ExtrusionMethod) method;
1341     int nameSize = 0;
1342     if ( load >> nameSize && nameSize > 0 )
1343     {
1344       _groupName.resize( nameSize );
1345       load.get( _groupName[0] ); // remove a white-space
1346       load.getline( &_groupName[0], nameSize + 1 );
1347     }
1348   }
1349   else {
1350     _isToIgnoreShapes = true; // old behavior
1351   }
1352   return load;
1353 } // --------------------------------------------------------------------------------
1354 bool StdMeshers_ViscousLayers::SetParametersByMesh(const SMESH_Mesh*   theMesh,
1355                                                    const TopoDS_Shape& theShape)
1356 {
1357   // TODO
1358   return false;
1359 } // --------------------------------------------------------------------------------
1360 SMESH_ComputeErrorPtr
1361 StdMeshers_ViscousLayers::CheckHypothesis(SMESH_Mesh&                          theMesh,
1362                                           const TopoDS_Shape&                  theShape,
1363                                           SMESH_Hypothesis::Hypothesis_Status& theStatus)
1364 {
1365   VISCOUS_3D::_ViscousBuilder builder;
1366   SMESH_ComputeErrorPtr err = builder.CheckHypotheses( theMesh, theShape );
1367   if ( err && !err->IsOK() )
1368     theStatus = SMESH_Hypothesis::HYP_INCOMPAT_HYPS;
1369   else
1370     theStatus = SMESH_Hypothesis::HYP_OK;
1371
1372   return err;
1373 }
1374 // --------------------------------------------------------------------------------
1375 bool StdMeshers_ViscousLayers::IsShapeWithLayers(int shapeIndex) const
1376 {
1377   bool isIn =
1378     ( std::find( _shapeIds.begin(), _shapeIds.end(), shapeIndex ) != _shapeIds.end() );
1379   return IsToIgnoreShapes() ? !isIn : isIn;
1380 }
1381
1382 // --------------------------------------------------------------------------------
1383 SMDS_MeshGroup* StdMeshers_ViscousLayers::CreateGroup( const std::string&  theName,
1384                                                        SMESH_Mesh&         theMesh,
1385                                                        SMDSAbs_ElementType theType)
1386 {
1387   SMESH_Group*      group = 0;
1388   SMDS_MeshGroup* groupDS = 0;
1389
1390   if ( theName.empty() )
1391     return groupDS;
1392        
1393   if ( SMESH_Mesh::GroupIteratorPtr grIt = theMesh.GetGroups() )
1394     while( grIt->more() && !group )
1395     {
1396       group = grIt->next();
1397       if ( !group ||
1398            group->GetGroupDS()->GetType() != theType ||
1399            group->GetName()               != theName ||
1400            !dynamic_cast< SMESHDS_Group* >( group->GetGroupDS() ))
1401         group = 0;
1402     }
1403   if ( !group )
1404     group = theMesh.AddGroup( theType, theName.c_str() );
1405
1406   groupDS = & dynamic_cast< SMESHDS_Group* >( group->GetGroupDS() )->SMDSGroup();
1407
1408   return groupDS;
1409 }
1410
1411 // END StdMeshers_ViscousLayers hypothesis
1412 //================================================================================
1413
1414 namespace VISCOUS_3D
1415 {
1416   gp_XYZ getEdgeDir( const TopoDS_Edge& E, const TopoDS_Vertex& fromV )
1417   {
1418     gp_Vec dir;
1419     double f,l;
1420     Handle(Geom_Curve) c = BRep_Tool::Curve( E, f, l );
1421     if ( c.IsNull() ) return gp_XYZ( Precision::Infinite(), 1e100, 1e100 );
1422     gp_Pnt p = BRep_Tool::Pnt( fromV );
1423     double distF = p.SquareDistance( c->Value( f ));
1424     double distL = p.SquareDistance( c->Value( l ));
1425     c->D1(( distF < distL ? f : l), p, dir );
1426     if ( distL < distF ) dir.Reverse();
1427     return dir.XYZ();
1428   }
1429   //--------------------------------------------------------------------------------
1430   gp_XYZ getEdgeDir( const TopoDS_Edge& E, const SMDS_MeshNode* atNode,
1431                      SMESH_MesherHelper& helper)
1432   {
1433     gp_Vec dir;
1434     double f,l; gp_Pnt p;
1435     Handle(Geom_Curve) c = BRep_Tool::Curve( E, f, l );
1436     if ( c.IsNull() ) return gp_XYZ( Precision::Infinite(), 1e100, 1e100 );
1437     double u = helper.GetNodeU( E, atNode );
1438     c->D1( u, p, dir );
1439     return dir.XYZ();
1440   }
1441   //--------------------------------------------------------------------------------
1442   gp_XYZ getFaceDir( const TopoDS_Face& F, const TopoDS_Vertex& fromV,
1443                      const SMDS_MeshNode* node, SMESH_MesherHelper& helper, bool& ok,
1444                      double* cosin=0);
1445   //--------------------------------------------------------------------------------
1446   gp_XYZ getFaceDir( const TopoDS_Face& F, const TopoDS_Edge& fromE,
1447                      const SMDS_MeshNode* node, SMESH_MesherHelper& helper, bool& ok)
1448   {
1449     double f,l;
1450     Handle(Geom_Curve) c = BRep_Tool::Curve( fromE, f, l );
1451     if ( c.IsNull() )
1452     {
1453       TopoDS_Vertex v = helper.IthVertex( 0, fromE );
1454       return getFaceDir( F, v, node, helper, ok );
1455     }
1456     gp_XY uv = helper.GetNodeUV( F, node, 0, &ok );
1457     Handle(Geom_Surface) surface = BRep_Tool::Surface( F );
1458     gp_Pnt p; gp_Vec du, dv, norm;
1459     surface->D1( uv.X(),uv.Y(), p, du,dv );
1460     norm = du ^ dv;
1461
1462     double u = helper.GetNodeU( fromE, node, 0, &ok );
1463     c->D1( u, p, du );
1464     TopAbs_Orientation o = helper.GetSubShapeOri( F.Oriented(TopAbs_FORWARD), fromE);
1465     if ( o == TopAbs_REVERSED )
1466       du.Reverse();
1467
1468     gp_Vec dir = norm ^ du;
1469
1470     if ( node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX &&
1471          helper.IsClosedEdge( fromE ))
1472     {
1473       if ( fabs(u-f) < fabs(u-l)) c->D1( l, p, dv );
1474       else                        c->D1( f, p, dv );
1475       if ( o == TopAbs_REVERSED )
1476         dv.Reverse();
1477       gp_Vec dir2 = norm ^ dv;
1478       dir = dir.Normalized() + dir2.Normalized();
1479     }
1480     return dir.XYZ();
1481   }
1482   //--------------------------------------------------------------------------------
1483   gp_XYZ getFaceDir( const TopoDS_Face& F, const TopoDS_Vertex& fromV,
1484                      const SMDS_MeshNode* node, SMESH_MesherHelper& helper,
1485                      bool& ok, double* cosin)
1486   {
1487     TopoDS_Face faceFrw = F;
1488     faceFrw.Orientation( TopAbs_FORWARD );
1489     //double f,l; TopLoc_Location loc;
1490     TopoDS_Edge edges[2]; // sharing a vertex
1491     size_t nbEdges = 0;
1492     {
1493       TopoDS_Vertex VV[2];
1494       TopExp_Explorer exp( faceFrw, TopAbs_EDGE );
1495       for ( ; exp.More() && nbEdges < 2; exp.Next() )
1496       {
1497         const TopoDS_Edge& e = TopoDS::Edge( exp.Current() );
1498         if ( SMESH_Algo::isDegenerated( e )) continue;
1499         TopExp::Vertices( e, VV[0], VV[1], /*CumOri=*/true );
1500         if ( VV[1].IsSame( fromV )) {
1501           nbEdges += edges[ 0 ].IsNull();
1502           edges[ 0 ] = e;
1503         }
1504         else if ( VV[0].IsSame( fromV )) {
1505           nbEdges += edges[ 1 ].IsNull();
1506           edges[ 1 ] = e;
1507         }
1508       }
1509     }
1510     gp_XYZ dir(0,0,0), edgeDir[2];
1511     if ( nbEdges == 2 )
1512     {
1513       // get dirs of edges going fromV
1514       ok = true;
1515       for ( size_t i = 0; i < nbEdges && ok; ++i )
1516       {
1517         edgeDir[i] = getEdgeDir( edges[i], fromV );
1518         double size2 = edgeDir[i].SquareModulus();
1519         if (( ok = size2 > numeric_limits<double>::min() ))
1520           edgeDir[i] /= sqrt( size2 );
1521       }
1522       if ( !ok ) return dir;
1523
1524       // get angle between the 2 edges
1525       gp_Vec faceNormal;
1526       double angle = helper.GetAngle( edges[0], edges[1], faceFrw, fromV, &faceNormal );
1527       if ( Abs( angle ) < 5 * M_PI/180 )
1528       {
1529         dir = ( faceNormal.XYZ() ^ edgeDir[0].Reversed()) + ( faceNormal.XYZ() ^ edgeDir[1] );
1530       }
1531       else
1532       {
1533         dir = edgeDir[0] + edgeDir[1];
1534         if ( angle < 0 )
1535           dir.Reverse();
1536       }
1537       if ( cosin ) {
1538         double angle = gp_Vec( edgeDir[0] ).Angle( dir );
1539         *cosin = Cos( angle );
1540       }
1541     }
1542     else if ( nbEdges == 1 )
1543     {
1544       dir = getFaceDir( faceFrw, edges[ edges[0].IsNull() ], node, helper, ok );
1545       if ( cosin ) *cosin = 1.;
1546     }
1547     else
1548     {
1549       ok = false;
1550     }
1551
1552     return dir;
1553   }
1554
1555   //================================================================================
1556   /*!
1557    * \brief Finds concave VERTEXes of a FACE
1558    */
1559   //================================================================================
1560
1561   bool getConcaveVertices( const TopoDS_Face&  F,
1562                            SMESH_MesherHelper& helper,
1563                            set< TGeomID >*     vertices = 0)
1564   {
1565     // check angles at VERTEXes
1566     TError error;
1567     TSideVector wires = StdMeshers_FaceSide::GetFaceWires( F, *helper.GetMesh(), 0, error );
1568     for ( size_t iW = 0; iW < wires.size(); ++iW )
1569     {
1570       const int nbEdges = wires[iW]->NbEdges();
1571       if ( nbEdges < 2 && SMESH_Algo::isDegenerated( wires[iW]->Edge(0)))
1572         continue;
1573       for ( int iE1 = 0; iE1 < nbEdges; ++iE1 )
1574       {
1575         if ( SMESH_Algo::isDegenerated( wires[iW]->Edge( iE1 ))) continue;
1576         int iE2 = ( iE1 + 1 ) % nbEdges;
1577         while ( SMESH_Algo::isDegenerated( wires[iW]->Edge( iE2 )))
1578           iE2 = ( iE2 + 1 ) % nbEdges;
1579         TopoDS_Vertex V = wires[iW]->FirstVertex( iE2 );
1580         double angle = helper.GetAngle( wires[iW]->Edge( iE1 ),
1581                                         wires[iW]->Edge( iE2 ), F, V );
1582         if ( angle < -5. * M_PI / 180. )
1583         {
1584           if ( !vertices )
1585             return true;
1586           vertices->insert( helper.GetMeshDS()->ShapeToIndex( V ));
1587         }
1588       }
1589     }
1590     return vertices ? !vertices->empty() : false;
1591   }
1592
1593   //================================================================================
1594   /*!
1595    * \brief Returns true if a FACE is bound by a concave EDGE
1596    */
1597   //================================================================================
1598
1599   bool isConcave( const TopoDS_Face&  F,
1600                   SMESH_MesherHelper& helper,
1601                   set< TGeomID >*     vertices = 0 )
1602   {
1603     bool isConcv = false;
1604     // if ( helper.Count( F, TopAbs_WIRE, /*useMap=*/false) > 1 )
1605     //   return true;
1606     gp_Vec2d drv1, drv2;
1607     gp_Pnt2d p;
1608     TopExp_Explorer eExp( F.Oriented( TopAbs_FORWARD ), TopAbs_EDGE );
1609     for ( ; eExp.More(); eExp.Next() )
1610     {
1611       const TopoDS_Edge& E = TopoDS::Edge( eExp.Current() );
1612       if ( SMESH_Algo::isDegenerated( E )) continue;
1613       // check if 2D curve is concave
1614       BRepAdaptor_Curve2d curve( E, F );
1615       const int nbIntervals = curve.NbIntervals( GeomAbs_C2 );
1616       TColStd_Array1OfReal intervals(1, nbIntervals + 1 );
1617       curve.Intervals( intervals, GeomAbs_C2 );
1618       bool isConvex = true;
1619       for ( int i = 1; i <= nbIntervals && isConvex; ++i )
1620       {
1621         double u1 = intervals( i );
1622         double u2 = intervals( i+1 );
1623         curve.D2( 0.5*( u1+u2 ), p, drv1, drv2 );
1624         double cross = drv1 ^ drv2;
1625         if ( E.Orientation() == TopAbs_REVERSED )
1626           cross = -cross;
1627         isConvex = ( cross > -1e-9 ); // 0.1 );
1628       }
1629       if ( !isConvex )
1630       {
1631         //cout << "Concave FACE " << helper.GetMeshDS()->ShapeToIndex( F ) << endl;
1632         isConcv = true;
1633         if ( vertices )
1634           break;
1635         else
1636           return true;
1637       }
1638     }
1639
1640     // check angles at VERTEXes
1641     if ( getConcaveVertices( F, helper, vertices ))
1642       isConcv = true;
1643
1644     return isConcv;
1645   }
1646
1647   //================================================================================
1648   /*!
1649    * \brief Computes minimal distance of face in-FACE nodes from an EDGE
1650    *  \param [in] face - the mesh face to treat
1651    *  \param [in] nodeOnEdge - a node on the EDGE
1652    *  \param [out] faceSize - the computed distance
1653    *  \return bool - true if faceSize computed
1654    */
1655   //================================================================================
1656
1657   bool getDistFromEdge( const SMDS_MeshElement* face,
1658                         const SMDS_MeshNode*    nodeOnEdge,
1659                         double &                faceSize )
1660   {
1661     faceSize = Precision::Infinite();
1662     bool done = false;
1663
1664     int nbN  = face->NbCornerNodes();
1665     int iOnE = face->GetNodeIndex( nodeOnEdge );
1666     int iNext[2] = { SMESH_MesherHelper::WrapIndex( iOnE+1, nbN ),
1667                      SMESH_MesherHelper::WrapIndex( iOnE-1, nbN ) };
1668     const SMDS_MeshNode* nNext[2] = { face->GetNode( iNext[0] ),
1669                                       face->GetNode( iNext[1] ) };
1670     gp_XYZ segVec, segEnd = SMESH_TNodeXYZ( nodeOnEdge ); // segment on EDGE
1671     double segLen = -1.;
1672     // look for two neighbor not in-FACE nodes of face
1673     for ( int i = 0; i < 2; ++i )
1674     {
1675       if (( nNext[i]->GetPosition()->GetDim() != 2 ) &&
1676           ( nodeOnEdge->GetPosition()->GetDim() == 0 || nNext[i]->GetID() < nodeOnEdge->GetID() ))
1677       {
1678         // look for an in-FACE node
1679         for ( int iN = 0; iN < nbN; ++iN )
1680         {
1681           if ( iN == iOnE || iN == iNext[i] )
1682             continue;
1683           SMESH_TNodeXYZ pInFace = face->GetNode( iN );
1684           gp_XYZ v = pInFace - segEnd;
1685           if ( segLen < 0 )
1686           {
1687             segVec = SMESH_TNodeXYZ( nNext[i] ) - segEnd;
1688             segLen = segVec.Modulus();
1689           }
1690           double distToSeg = v.Crossed( segVec ).Modulus() / segLen;
1691           faceSize = Min( faceSize, distToSeg );
1692           done = true;
1693         }
1694         segLen = -1;
1695       }
1696     }
1697     return done;
1698   }
1699   //================================================================================
1700   /*!
1701    * \brief Return direction of axis or revolution of a surface
1702    */
1703   //================================================================================
1704
1705   bool getRovolutionAxis( const Adaptor3d_Surface& surface,
1706                           gp_Dir &                 axis )
1707   {
1708     switch ( surface.GetType() ) {
1709     case GeomAbs_Cone:
1710     {
1711       gp_Cone cone = surface.Cone();
1712       axis = cone.Axis().Direction();
1713       break;
1714     }
1715     case GeomAbs_Sphere:
1716     {
1717       gp_Sphere sphere = surface.Sphere();
1718       axis = sphere.Position().Direction();
1719       break;
1720     }
1721     case GeomAbs_SurfaceOfRevolution:
1722     {
1723       axis = surface.AxeOfRevolution().Direction();
1724       break;
1725     }
1726     //case GeomAbs_SurfaceOfExtrusion:
1727     case GeomAbs_OffsetSurface:
1728     {
1729       Handle(Adaptor3d_HSurface) base = surface.BasisSurface();
1730       return getRovolutionAxis( base->Surface(), axis );
1731     }
1732     default: return false;
1733     }
1734     return true;
1735   }
1736
1737   //--------------------------------------------------------------------------------
1738   // DEBUG. Dump intermediate node positions into a python script
1739   // HOWTO use: run python commands written in a console to see
1740   //  construction steps of viscous layers
1741 #ifdef __myDEBUG
1742   ostream* py;
1743   int      theNbPyFunc;
1744   struct PyDump
1745   {
1746     PyDump(SMESH_Mesh& m) {
1747       int tag = 3 + m.GetId();
1748       const char* fname = "/tmp/viscous.py";
1749       cout << "execfile('"<<fname<<"')"<<endl;
1750       py = _pyStream = new ofstream(fname);
1751       *py << "import SMESH" << endl
1752           << "from salome.smesh import smeshBuilder" << endl
1753           << "smesh  = smeshBuilder.New()" << endl
1754           << "meshSO = salome.myStudy.FindObjectID('0:1:2:" << tag <<"')" << endl
1755           << "mesh   = smesh.Mesh( meshSO.GetObject() )"<<endl;
1756       theNbPyFunc = 0;
1757     }
1758     void Finish() {
1759       if (py) {
1760         *py << "mesh.GroupOnFilter(SMESH.VOLUME,'Viscous Prisms',"
1761           "smesh.GetFilter(SMESH.VOLUME,SMESH.FT_ElemGeomType,'=',SMESH.Geom_PENTA))"<<endl;
1762         *py << "mesh.GroupOnFilter(SMESH.VOLUME,'Neg Volumes',"
1763           "smesh.GetFilter(SMESH.VOLUME,SMESH.FT_Volume3D,'<',0))"<<endl;
1764       }
1765       delete py; py=0;
1766     }
1767     ~PyDump() { Finish(); cout << "NB FUNCTIONS: " << theNbPyFunc << endl; }
1768     struct MyStream : public ostream
1769     {
1770       template <class T> ostream & operator<<( const T &anything ) { return *this ; }
1771     };
1772     void Pause() { py = &_mystream; }
1773     void Resume() { py = _pyStream; }
1774     MyStream _mystream;
1775     ostream* _pyStream;
1776   };
1777 #define dumpFunction(f) { _dumpFunction(f, __LINE__);}
1778 #define dumpMove(n)     { _dumpMove(n, __LINE__);}
1779 #define dumpMoveComm(n,txt) { _dumpMove(n, __LINE__, txt);}
1780 #define dumpCmd(txt)    { _dumpCmd(txt, __LINE__);}
1781   void _dumpFunction(const string& fun, int ln)
1782   { if (py) *py<< "def "<<fun<<"(): # "<< ln <<endl; cout<<fun<<"()"<<endl; ++theNbPyFunc; }
1783   void _dumpMove(const SMDS_MeshNode* n, int ln, const char* txt="")
1784   { if (py) *py<< "  mesh.MoveNode( "<<n->GetID()<< ", "<< n->X()
1785                << ", "<<n->Y()<<", "<< n->Z()<< ")\t\t # "<< ln <<" "<< txt << endl; }
1786   void _dumpCmd(const string& txt, int ln)
1787   { if (py) *py<< "  "<<txt<<" # "<< ln <<endl; }
1788   void dumpFunctionEnd()
1789   { if (py) *py<< "  return"<< endl; }
1790   void dumpChangeNodes( const SMDS_MeshElement* f )
1791   { if (py) { *py<< "  mesh.ChangeElemNodes( " << f->GetID()<<", [";
1792       for ( int i=1; i < f->NbNodes(); ++i ) *py << f->GetNode(i-1)->GetID()<<", ";
1793       *py << f->GetNode( f->NbNodes()-1 )->GetID() << " ])"<< endl; }}
1794 #define debugMsg( txt ) { cout << "# "<< txt << " (line: " << __LINE__ << ")" << endl; }
1795
1796 #else
1797
1798   struct PyDump { PyDump(SMESH_Mesh&) {} void Finish() {} void Pause() {} void Resume() {} };
1799 #define dumpFunction(f) f
1800 #define dumpMove(n)
1801 #define dumpMoveComm(n,txt)
1802 #define dumpCmd(txt)
1803 #define dumpFunctionEnd()
1804 #define dumpChangeNodes(f) { if(f) {} } // prevent "unused variable 'f'" warning
1805 #define debugMsg( txt ) {}
1806
1807 #endif
1808 }
1809
1810 using namespace VISCOUS_3D;
1811
1812 //================================================================================
1813 /*!
1814  * \brief Constructor of _ViscousBuilder
1815  */
1816 //================================================================================
1817
1818 _ViscousBuilder::_ViscousBuilder()
1819 {
1820   _error = SMESH_ComputeError::New(COMPERR_OK);
1821   _tmpFaceID = 0;
1822 }
1823
1824 //================================================================================
1825 /*!
1826  * \brief Stores error description and returns false
1827  */
1828 //================================================================================
1829
1830 bool _ViscousBuilder::error(const string& text, int solidId )
1831 {
1832   const string prefix = string("Viscous layers builder: ");
1833   _error->myName    = COMPERR_ALGO_FAILED;
1834   _error->myComment = prefix + text;
1835   if ( _mesh )
1836   {
1837     SMESH_subMesh* sm = _mesh->GetSubMeshContaining( solidId );
1838     if ( !sm && !_sdVec.empty() )
1839       sm = _mesh->GetSubMeshContaining( solidId = _sdVec[0]._index );
1840     if ( sm && sm->GetSubShape().ShapeType() == TopAbs_SOLID )
1841     {
1842       SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1843       if ( smError && smError->myAlgo )
1844         _error->myAlgo = smError->myAlgo;
1845       smError = _error;
1846       sm->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
1847     }
1848     // set KO to all solids
1849     for ( size_t i = 0; i < _sdVec.size(); ++i )
1850     {
1851       if ( _sdVec[i]._index == solidId )
1852         continue;
1853       sm = _mesh->GetSubMesh( _sdVec[i]._solid );
1854       if ( !sm->IsEmpty() )
1855         continue;
1856       SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1857       if ( !smError || smError->IsOK() )
1858       {
1859         smError = SMESH_ComputeError::New( COMPERR_ALGO_FAILED, prefix + "failed");
1860         sm->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
1861       }
1862     }
1863   }
1864   makeGroupOfLE(); // debug
1865
1866   return false;
1867 }
1868
1869 //================================================================================
1870 /*!
1871  * \brief At study restoration, restore event listeners used to clear an inferior
1872  *  dim sub-mesh modified by viscous layers
1873  */
1874 //================================================================================
1875
1876 void _ViscousBuilder::RestoreListeners()
1877 {
1878   // TODO
1879 }
1880
1881 //================================================================================
1882 /*!
1883  * \brief computes SMESH_ProxyMesh::SubMesh::_n2n
1884  */
1885 //================================================================================
1886
1887 bool _ViscousBuilder::MakeN2NMap( _MeshOfSolid* pm )
1888 {
1889   SMESH_subMesh* solidSM = pm->mySubMeshes.front();
1890   TopExp_Explorer fExp( solidSM->GetSubShape(), TopAbs_FACE );
1891   for ( ; fExp.More(); fExp.Next() )
1892   {
1893     SMESHDS_SubMesh* srcSmDS = pm->GetMeshDS()->MeshElements( fExp.Current() );
1894     const SMESH_ProxyMesh::SubMesh* prxSmDS = pm->GetProxySubMesh( fExp.Current() );
1895
1896     if ( !srcSmDS || !prxSmDS || !srcSmDS->NbElements() || !prxSmDS->NbElements() )
1897       continue;
1898     if ( srcSmDS->GetElements()->next() == prxSmDS->GetElements()->next())
1899       continue;
1900
1901     if ( srcSmDS->NbElements() != prxSmDS->NbElements() )
1902       return error( "Different nb elements in a source and a proxy sub-mesh", solidSM->GetId());
1903
1904     SMDS_ElemIteratorPtr srcIt = srcSmDS->GetElements();
1905     SMDS_ElemIteratorPtr prxIt = prxSmDS->GetElements();
1906     while( prxIt->more() )
1907     {
1908       const SMDS_MeshElement* fSrc = srcIt->next();
1909       const SMDS_MeshElement* fPrx = prxIt->next();
1910       if ( fSrc->NbNodes() != fPrx->NbNodes())
1911         return error( "Different elements in a source and a proxy sub-mesh", solidSM->GetId());
1912       for ( int i = 0 ; i < fPrx->NbNodes(); ++i )
1913         pm->setNode2Node( fSrc->GetNode(i), fPrx->GetNode(i), prxSmDS );
1914     }
1915   }
1916   pm->_n2nMapComputed = true;
1917   return true;
1918 }
1919
1920 //================================================================================
1921 /*!
1922  * \brief Does its job
1923  */
1924 //================================================================================
1925
1926 SMESH_ComputeErrorPtr _ViscousBuilder::Compute(SMESH_Mesh&         theMesh,
1927                                                const TopoDS_Shape& theShape)
1928 {
1929   _mesh = & theMesh;
1930
1931   // check if proxy mesh already computed
1932   TopExp_Explorer exp( theShape, TopAbs_SOLID );
1933   if ( !exp.More() )
1934     return error("No SOLID's in theShape"), _error;
1935
1936   if ( _ViscousListener::GetSolidMesh( _mesh, exp.Current(), /*toCreate=*/false))
1937     return SMESH_ComputeErrorPtr(); // everything already computed
1938
1939   PyDump debugDump( theMesh );
1940   _pyDump = &debugDump;
1941
1942   // TODO: ignore already computed SOLIDs 
1943   if ( !findSolidsWithLayers())
1944     return _error;
1945
1946   if ( !findFacesWithLayers() )
1947     return _error;
1948
1949   for ( size_t i = 0; i < _sdVec.size(); ++i )
1950   {
1951     size_t iSD = 0;
1952     for ( iSD = 0; iSD < _sdVec.size(); ++iSD ) // find next SOLID to compute
1953       if ( _sdVec[iSD]._before.IsEmpty() &&
1954            !_sdVec[iSD]._solid.IsNull() &&
1955            _sdVec[iSD]._n2eMap.empty() )
1956         break;
1957
1958     if ( ! makeLayer(_sdVec[iSD]) )   // create _LayerEdge's
1959       return _error;
1960
1961     if ( _sdVec[iSD]._n2eMap.size() == 0 ) // no layers in a SOLID
1962     {
1963       _sdVec[iSD]._solid.Nullify();
1964       continue;
1965     }
1966
1967     if ( ! inflate(_sdVec[iSD]) )     // increase length of _LayerEdge's
1968       return _error;
1969
1970     if ( ! refine(_sdVec[iSD]) )      // create nodes and prisms
1971       return _error;
1972
1973     if ( ! shrink(_sdVec[iSD]) )      // shrink 2D mesh on FACEs w/o layer
1974       return _error;
1975
1976     addBoundaryElements(_sdVec[iSD]); // create quadrangles on prism bare sides
1977
1978     const TopoDS_Shape& solid = _sdVec[iSD]._solid;
1979     for ( iSD = 0; iSD < _sdVec.size(); ++iSD )
1980       _sdVec[iSD]._before.Remove( solid );
1981   }
1982
1983   makeGroupOfLE(); // debug
1984   debugDump.Finish();
1985
1986   return _error;
1987 }
1988
1989 //================================================================================
1990 /*!
1991  * \brief Check validity of hypotheses
1992  */
1993 //================================================================================
1994
1995 SMESH_ComputeErrorPtr _ViscousBuilder::CheckHypotheses( SMESH_Mesh&         mesh,
1996                                                         const TopoDS_Shape& shape )
1997 {
1998   _mesh = & mesh;
1999
2000   if ( _ViscousListener::GetSolidMesh( _mesh, shape, /*toCreate=*/false))
2001     return SMESH_ComputeErrorPtr(); // everything already computed
2002
2003
2004   findSolidsWithLayers();
2005   bool ok = findFacesWithLayers( true );
2006
2007   // remove _MeshOfSolid's of _SolidData's
2008   for ( size_t i = 0; i < _sdVec.size(); ++i )
2009     _ViscousListener::RemoveSolidMesh( _mesh, _sdVec[i]._solid );
2010
2011   if ( !ok )
2012     return _error;
2013
2014   return SMESH_ComputeErrorPtr();
2015 }
2016
2017 //================================================================================
2018 /*!
2019  * \brief Finds SOLIDs to compute using viscous layers. Fills _sdVec
2020  */
2021 //================================================================================
2022
2023 bool _ViscousBuilder::findSolidsWithLayers()
2024 {
2025   // get all solids
2026   TopTools_IndexedMapOfShape allSolids;
2027   TopExp::MapShapes( _mesh->GetShapeToMesh(), TopAbs_SOLID, allSolids );
2028   _sdVec.reserve( allSolids.Extent());
2029
2030   SMESH_HypoFilter filter;
2031   for ( int i = 1; i <= allSolids.Extent(); ++i )
2032   {
2033     // find StdMeshers_ViscousLayers hyp assigned to the i-th solid
2034     SMESH_subMesh* sm = _mesh->GetSubMesh( allSolids(i) );
2035     if ( sm->GetSubMeshDS() && sm->GetSubMeshDS()->NbElements() > 0 )
2036       continue; // solid is already meshed
2037     SMESH_Algo* algo = sm->GetAlgo();
2038     if ( !algo ) continue;
2039     // TODO: check if algo is hidden
2040     const list <const SMESHDS_Hypothesis *> & allHyps =
2041       algo->GetUsedHypothesis(*_mesh, allSolids(i), /*ignoreAuxiliary=*/false);
2042     _SolidData* soData = 0;
2043     list< const SMESHDS_Hypothesis *>::const_iterator hyp = allHyps.begin();
2044     const StdMeshers_ViscousLayers* viscHyp = 0;
2045     for ( ; hyp != allHyps.end(); ++hyp )
2046       if (( viscHyp = dynamic_cast<const StdMeshers_ViscousLayers*>( *hyp )))
2047       {
2048         TopoDS_Shape hypShape;
2049         filter.Init( filter.Is( viscHyp ));
2050         _mesh->GetHypothesis( allSolids(i), filter, true, &hypShape );
2051
2052         if ( !soData )
2053         {
2054           _MeshOfSolid* proxyMesh = _ViscousListener::GetSolidMesh( _mesh,
2055                                                                     allSolids(i),
2056                                                                     /*toCreate=*/true);
2057           _sdVec.push_back( _SolidData( allSolids(i), proxyMesh ));
2058           soData = & _sdVec.back();
2059           soData->_index = getMeshDS()->ShapeToIndex( allSolids(i));
2060           soData->_helper = new SMESH_MesherHelper( *_mesh );
2061           soData->_helper->SetSubShape( allSolids(i) );
2062           _solids.Add( allSolids(i) );
2063         }
2064         soData->_hyps.push_back( viscHyp );
2065         soData->_hypShapes.push_back( hypShape );
2066       }
2067   }
2068   if ( _sdVec.empty() )
2069     return error
2070       ( SMESH_Comment(StdMeshers_ViscousLayers::GetHypType()) << " hypothesis not found",0);
2071
2072   return true;
2073 }
2074
2075 //================================================================================
2076 /*!
2077  * \brief Set a _SolidData to be computed before another
2078  */
2079 //================================================================================
2080
2081 bool _ViscousBuilder::setBefore( _SolidData& solidBefore, _SolidData& solidAfter )
2082 {
2083   // check possibility to set this order; get all solids before solidBefore
2084   TopTools_IndexedMapOfShape allSolidsBefore;
2085   allSolidsBefore.Add( solidBefore._solid );
2086   for ( int i = 1; i <= allSolidsBefore.Extent(); ++i )
2087   {
2088     int iSD = _solids.FindIndex( allSolidsBefore(i) );
2089     if ( iSD )
2090     {
2091       TopTools_MapIteratorOfMapOfShape soIt( _sdVec[ iSD-1 ]._before );
2092       for ( ; soIt.More(); soIt.Next() )
2093         allSolidsBefore.Add( soIt.Value() );
2094     }
2095   }
2096   if ( allSolidsBefore.Contains( solidAfter._solid ))
2097     return false;
2098
2099   for ( int i = 1; i <= allSolidsBefore.Extent(); ++i )
2100     solidAfter._before.Add( allSolidsBefore(i) );
2101
2102   return true;
2103 }
2104
2105 //================================================================================
2106 /*!
2107  * \brief
2108  */
2109 //================================================================================
2110
2111 bool _ViscousBuilder::findFacesWithLayers(const bool onlyWith)
2112 {
2113   SMESH_MesherHelper helper( *_mesh );
2114   TopExp_Explorer exp;
2115
2116   // collect all faces-to-ignore defined by hyp
2117   for ( size_t i = 0; i < _sdVec.size(); ++i )
2118   {
2119     // get faces-to-ignore defined by each hyp
2120     typedef const StdMeshers_ViscousLayers* THyp;
2121     typedef std::pair< set<TGeomID>, THyp > TFacesOfHyp;
2122     list< TFacesOfHyp > ignoreFacesOfHyps;
2123     list< THyp >::iterator              hyp = _sdVec[i]._hyps.begin();
2124     list< TopoDS_Shape >::iterator hypShape = _sdVec[i]._hypShapes.begin();
2125     for ( ; hyp != _sdVec[i]._hyps.end(); ++hyp, ++hypShape )
2126     {
2127       ignoreFacesOfHyps.push_back( TFacesOfHyp( set<TGeomID>(), *hyp ));
2128       getIgnoreFaces( _sdVec[i]._solid, *hyp, *hypShape, ignoreFacesOfHyps.back().first );
2129     }
2130
2131     // fill _SolidData::_face2hyp and check compatibility of hypotheses
2132     const int nbHyps = _sdVec[i]._hyps.size();
2133     if ( nbHyps > 1 )
2134     {
2135       // check if two hypotheses define different parameters for the same FACE
2136       list< TFacesOfHyp >::iterator igFacesOfHyp;
2137       for ( exp.Init( _sdVec[i]._solid, TopAbs_FACE ); exp.More(); exp.Next() )
2138       {
2139         const TGeomID faceID = getMeshDS()->ShapeToIndex( exp.Current() );
2140         THyp hyp = 0;
2141         igFacesOfHyp = ignoreFacesOfHyps.begin();
2142         for ( ; igFacesOfHyp != ignoreFacesOfHyps.end(); ++igFacesOfHyp )
2143           if ( ! igFacesOfHyp->first.count( faceID ))
2144           {
2145             if ( hyp )
2146               return error(SMESH_Comment("Several hypotheses define "
2147                                          "Viscous Layers on the face #") << faceID );
2148             hyp = igFacesOfHyp->second;
2149           }
2150         if ( hyp )
2151           _sdVec[i]._face2hyp.insert( make_pair( faceID, hyp ));
2152         else
2153           _sdVec[i]._ignoreFaceIds.insert( faceID );
2154       }
2155
2156       // check if two hypotheses define different number of viscous layers for
2157       // adjacent faces of a solid
2158       set< int > nbLayersSet;
2159       igFacesOfHyp = ignoreFacesOfHyps.begin();
2160       for ( ; igFacesOfHyp != ignoreFacesOfHyps.end(); ++igFacesOfHyp )
2161       {
2162         nbLayersSet.insert( igFacesOfHyp->second->GetNumberLayers() );
2163       }
2164       if ( nbLayersSet.size() > 1 )
2165       {
2166         for ( exp.Init( _sdVec[i]._solid, TopAbs_EDGE ); exp.More(); exp.Next() )
2167         {
2168           PShapeIteratorPtr fIt = helper.GetAncestors( exp.Current(), *_mesh, TopAbs_FACE );
2169           THyp hyp1 = 0, hyp2 = 0;
2170           while( const TopoDS_Shape* face = fIt->next() )
2171           {
2172             const TGeomID faceID = getMeshDS()->ShapeToIndex( *face );
2173             map< TGeomID, THyp >::iterator f2h = _sdVec[i]._face2hyp.find( faceID );
2174             if ( f2h != _sdVec[i]._face2hyp.end() )
2175             {
2176               ( hyp1 ? hyp2 : hyp1 ) = f2h->second;
2177             }
2178           }
2179           if ( hyp1 && hyp2 &&
2180                hyp1->GetNumberLayers() != hyp2->GetNumberLayers() )
2181           {
2182             return error("Two hypotheses define different number of "
2183                          "viscous layers on adjacent faces");
2184           }
2185         }
2186       }
2187     } // if ( nbHyps > 1 )
2188     else
2189     {
2190       _sdVec[i]._ignoreFaceIds.swap( ignoreFacesOfHyps.back().first );
2191     }
2192   } // loop on _sdVec
2193
2194   if ( onlyWith ) // is called to check hypotheses compatibility only
2195     return true;
2196
2197   // fill _SolidData::_reversedFaceIds
2198   for ( size_t i = 0; i < _sdVec.size(); ++i )
2199   {
2200     exp.Init( _sdVec[i]._solid.Oriented( TopAbs_FORWARD ), TopAbs_FACE );
2201     for ( ; exp.More(); exp.Next() )
2202     {
2203       const TopoDS_Face& face = TopoDS::Face( exp.Current() );
2204       const TGeomID faceID = getMeshDS()->ShapeToIndex( face );
2205       if ( //!sdVec[i]._ignoreFaceIds.count( faceID ) &&
2206           helper.NbAncestors( face, *_mesh, TopAbs_SOLID ) > 1 &&
2207           helper.IsReversedSubMesh( face ))
2208       {
2209         _sdVec[i]._reversedFaceIds.insert( faceID );
2210       }
2211     }
2212   }
2213
2214   // Find FACEs to shrink mesh on (solution 2 in issue 0020832): fill in _shrinkShape2Shape
2215   TopTools_IndexedMapOfShape shapes;
2216   std::string structAlgoName = "Hexa_3D";
2217   for ( size_t i = 0; i < _sdVec.size(); ++i )
2218   {
2219     shapes.Clear();
2220     TopExp::MapShapes(_sdVec[i]._solid, TopAbs_EDGE, shapes);
2221     for ( int iE = 1; iE <= shapes.Extent(); ++iE )
2222     {
2223       const TopoDS_Shape& edge = shapes(iE);
2224       // find 2 FACEs sharing an EDGE
2225       TopoDS_Shape FF[2];
2226       PShapeIteratorPtr fIt = helper.GetAncestors(edge, *_mesh, TopAbs_FACE, &_sdVec[i]._solid);
2227       while ( fIt->more())
2228       {
2229         const TopoDS_Shape* f = fIt->next();
2230         FF[ int( !FF[0].IsNull()) ] = *f;
2231       }
2232       if( FF[1].IsNull() ) continue; // seam edge can be shared by 1 FACE only
2233
2234       // check presence of layers on them
2235       int ignore[2];
2236       for ( int j = 0; j < 2; ++j )
2237         ignore[j] = _sdVec[i]._ignoreFaceIds.count( getMeshDS()->ShapeToIndex( FF[j] ));
2238       if ( ignore[0] == ignore[1] )
2239         continue; // nothing interesting
2240       TopoDS_Shape fWOL = FF[ ignore[0] ? 0 : 1 ];
2241
2242       // add EDGE to maps
2243       if ( !fWOL.IsNull())
2244       {
2245         TGeomID edgeInd = getMeshDS()->ShapeToIndex( edge );
2246         _sdVec[i]._shrinkShape2Shape.insert( make_pair( edgeInd, fWOL ));
2247       }
2248     }
2249   }
2250
2251   // Find the SHAPE along which to inflate _LayerEdge based on VERTEX
2252
2253   for ( size_t i = 0; i < _sdVec.size(); ++i )
2254   {
2255     shapes.Clear();
2256     TopExp::MapShapes(_sdVec[i]._solid, TopAbs_VERTEX, shapes);
2257     for ( int iV = 1; iV <= shapes.Extent(); ++iV )
2258     {
2259       const TopoDS_Shape& vertex = shapes(iV);
2260       // find faces WOL sharing the vertex
2261       vector< TopoDS_Shape > facesWOL;
2262       size_t totalNbFaces = 0;
2263       PShapeIteratorPtr fIt = helper.GetAncestors(vertex, *_mesh, TopAbs_FACE, &_sdVec[i]._solid );
2264       while ( fIt->more())
2265       {
2266         const TopoDS_Shape* f = fIt->next();
2267         totalNbFaces++;
2268         const int fID = getMeshDS()->ShapeToIndex( *f );
2269         if ( _sdVec[i]._ignoreFaceIds.count ( fID ) /*&& !_sdVec[i]._noShrinkShapes.count( fID )*/)
2270           facesWOL.push_back( *f );
2271       }
2272       if ( facesWOL.size() == totalNbFaces || facesWOL.empty() )
2273         continue; // no layers at this vertex or no WOL
2274       TGeomID vInd = getMeshDS()->ShapeToIndex( vertex );
2275       switch ( facesWOL.size() )
2276       {
2277       case 1:
2278       {
2279         helper.SetSubShape( facesWOL[0] );
2280         if ( helper.IsRealSeam( vInd )) // inflate along a seam edge?
2281         {
2282           TopoDS_Shape seamEdge;
2283           PShapeIteratorPtr eIt = helper.GetAncestors(vertex, *_mesh, TopAbs_EDGE);
2284           while ( eIt->more() && seamEdge.IsNull() )
2285           {
2286             const TopoDS_Shape* e = eIt->next();
2287             if ( helper.IsRealSeam( *e ) )
2288               seamEdge = *e;
2289           }
2290           if ( !seamEdge.IsNull() )
2291           {
2292             _sdVec[i]._shrinkShape2Shape.insert( make_pair( vInd, seamEdge ));
2293             break;
2294           }
2295         }
2296         _sdVec[i]._shrinkShape2Shape.insert( make_pair( vInd, facesWOL[0] ));
2297         break;
2298       }
2299       case 2:
2300       {
2301         // find an edge shared by 2 faces
2302         PShapeIteratorPtr eIt = helper.GetAncestors(vertex, *_mesh, TopAbs_EDGE);
2303         while ( eIt->more())
2304         {
2305           const TopoDS_Shape* e = eIt->next();
2306           if ( helper.IsSubShape( *e, facesWOL[0]) &&
2307                helper.IsSubShape( *e, facesWOL[1]))
2308           {
2309             _sdVec[i]._shrinkShape2Shape.insert( make_pair( vInd, *e )); break;
2310           }
2311         }
2312         break;
2313       }
2314       default:
2315         return error("Not yet supported case", _sdVec[i]._index);
2316       }
2317     }
2318   }
2319
2320   // Add to _noShrinkShapes sub-shapes of FACE's that can't be shrinked since
2321   // the algo of the SOLID sharing the FACE does not support it or for other reasons
2322   set< string > notSupportAlgos; notSupportAlgos.insert( structAlgoName );
2323   for ( size_t i = 0; i < _sdVec.size(); ++i )
2324   {
2325     map< TGeomID, TopoDS_Shape >::iterator e2f = _sdVec[i]._shrinkShape2Shape.begin();
2326     for ( ; e2f != _sdVec[i]._shrinkShape2Shape.end(); ++e2f )
2327     {
2328       const TopoDS_Shape& fWOL = e2f->second;
2329       const TGeomID     edgeID = e2f->first;
2330       TGeomID           faceID = getMeshDS()->ShapeToIndex( fWOL );
2331       TopoDS_Shape        edge = getMeshDS()->IndexToShape( edgeID );
2332       if ( edge.ShapeType() != TopAbs_EDGE )
2333         continue; // shrink shape is VERTEX
2334
2335       TopoDS_Shape solid;
2336       PShapeIteratorPtr soIt = helper.GetAncestors(fWOL, *_mesh, TopAbs_SOLID);
2337       while ( soIt->more() && solid.IsNull() )
2338       {
2339         const TopoDS_Shape* so = soIt->next();
2340         if ( !so->IsSame( _sdVec[i]._solid ))
2341           solid = *so;
2342       }
2343       if ( solid.IsNull() )
2344         continue;
2345
2346       bool noShrinkE = false;
2347       SMESH_Algo*  algo = _mesh->GetSubMesh( solid )->GetAlgo();
2348       bool isStructured = ( algo && algo->GetName() == structAlgoName );
2349       size_t     iSolid = _solids.FindIndex( solid ) - 1;
2350       if ( iSolid < _sdVec.size() && _sdVec[ iSolid ]._ignoreFaceIds.count( faceID ))
2351       {
2352         // the adjacent SOLID has NO layers on fWOL;
2353         // shrink allowed if
2354         // - there are layers on the EDGE in the adjacent SOLID
2355         // - there are NO layers in the adjacent SOLID && algo is unstructured and computed later
2356         bool hasWLAdj = (_sdVec[iSolid]._shrinkShape2Shape.count( edgeID ));
2357         bool shrinkAllowed = (( hasWLAdj ) ||
2358                               ( !isStructured && setBefore( _sdVec[ i ], _sdVec[ iSolid ] )));
2359         noShrinkE = !shrinkAllowed;
2360       }
2361       else if ( iSolid < _sdVec.size() )
2362       {
2363         // the adjacent SOLID has layers on fWOL;
2364         // check if SOLID's mesh is unstructured and then try to set it
2365         // to be computed after the i-th solid
2366         if ( isStructured || !setBefore( _sdVec[ i ], _sdVec[ iSolid ] ))
2367           noShrinkE = true; // don't shrink fWOL
2368       }
2369       else
2370       {
2371         // the adjacent SOLID has NO layers at all
2372         noShrinkE = isStructured;
2373       }
2374
2375       if ( noShrinkE )
2376       {
2377         _sdVec[i]._noShrinkShapes.insert( edgeID );
2378
2379         // check if there is a collision with to-shrink-from EDGEs in iSolid
2380         // if ( iSolid < _sdVec.size() )
2381         // {
2382         //   shapes.Clear();
2383         //   TopExp::MapShapes( fWOL, TopAbs_EDGE, shapes);
2384         //   for ( int iE = 1; iE <= shapes.Extent(); ++iE )
2385         //   {
2386         //     const TopoDS_Edge& E = TopoDS::Edge( shapes( iE ));
2387         //     const TGeomID    eID = getMeshDS()->ShapeToIndex( E );
2388         //     if ( eID == edgeID ||
2389         //          !_sdVec[iSolid]._shrinkShape2Shape.count( eID ) ||
2390         //          _sdVec[i]._noShrinkShapes.count( eID ))
2391         //       continue;
2392         //     for ( int is1st = 0; is1st < 2; ++is1st )
2393         //     {
2394         //       TopoDS_Vertex V = helper.IthVertex( is1st, E );
2395         //       if ( _sdVec[i]._noShrinkShapes.count( getMeshDS()->ShapeToIndex( V ) ))
2396         //       {
2397         //         return error("No way to make a conformal mesh with "
2398         //                      "the given set of faces with layers", _sdVec[i]._index);
2399         //       }
2400         //     }
2401         //   }
2402         // }
2403       }
2404
2405       // add VERTEXes of the edge in _noShrinkShapes, which is necessary if
2406       // _shrinkShape2Shape is different in the adjacent SOLID
2407       for ( TopoDS_Iterator vIt( edge ); vIt.More(); vIt.Next() )
2408       {
2409         TGeomID vID = getMeshDS()->ShapeToIndex( vIt.Value() );
2410         bool noShrinkV = false, noShrinkIfAdjMeshed = false;
2411
2412         if ( iSolid < _sdVec.size() )
2413         {
2414           if ( _sdVec[ iSolid ]._ignoreFaceIds.count( faceID ))
2415           {
2416             map< TGeomID, TopoDS_Shape >::iterator i2S, i2SAdj;
2417             i2S    = _sdVec[i     ]._shrinkShape2Shape.find( vID );
2418             i2SAdj = _sdVec[iSolid]._shrinkShape2Shape.find( vID );
2419             if ( i2SAdj == _sdVec[iSolid]._shrinkShape2Shape.end() )
2420               noShrinkV = (( isStructured ) ||
2421                            ( noShrinkIfAdjMeshed = i2S->second.ShapeType() == TopAbs_EDGE ));
2422             else
2423               noShrinkV = ( ! i2S->second.IsSame( i2SAdj->second ));
2424           }
2425           else
2426           {
2427             noShrinkV = noShrinkE;
2428           }
2429         }
2430         else
2431         {
2432           // the adjacent SOLID has NO layers at all
2433           if ( isStructured )
2434           {
2435             noShrinkV = true;
2436           }
2437           else
2438           {
2439             noShrinkV = noShrinkIfAdjMeshed =
2440               ( _sdVec[i]._shrinkShape2Shape[ vID ].ShapeType() == TopAbs_EDGE );
2441           }
2442         }
2443
2444         if ( noShrinkV && noShrinkIfAdjMeshed )
2445         {
2446           // noShrinkV if FACEs in the adjacent SOLID are meshed
2447           PShapeIteratorPtr fIt = helper.GetAncestors( _sdVec[i]._shrinkShape2Shape[ vID ],
2448                                                        *_mesh, TopAbs_FACE, &solid );
2449           while ( fIt->more() )
2450           {
2451             const TopoDS_Shape* f = fIt->next();
2452             if ( !f->IsSame( fWOL ))
2453             {
2454               noShrinkV = ! _mesh->GetSubMesh( *f )->IsEmpty();
2455               break;
2456             }
2457           }
2458         }
2459         if ( noShrinkV )
2460           _sdVec[i]._noShrinkShapes.insert( vID );
2461       }
2462
2463     } // loop on _sdVec[i]._shrinkShape2Shape
2464   } // loop on _sdVec to fill in _SolidData::_noShrinkShapes
2465
2466
2467     // add FACEs of other SOLIDs to _ignoreFaceIds
2468   for ( size_t i = 0; i < _sdVec.size(); ++i )
2469   {
2470     shapes.Clear();
2471     TopExp::MapShapes(_sdVec[i]._solid, TopAbs_FACE, shapes);
2472
2473     for ( exp.Init( _mesh->GetShapeToMesh(), TopAbs_FACE ); exp.More(); exp.Next() )
2474     {
2475       if ( !shapes.Contains( exp.Current() ))
2476         _sdVec[i]._ignoreFaceIds.insert( getMeshDS()->ShapeToIndex( exp.Current() ));
2477     }
2478   }
2479
2480   return true;
2481 }
2482
2483 //================================================================================
2484 /*!
2485  * \brief Finds FACEs w/o layers for a given SOLID by an hypothesis
2486  */
2487 //================================================================================
2488
2489 void _ViscousBuilder::getIgnoreFaces(const TopoDS_Shape&             solid,
2490                                      const StdMeshers_ViscousLayers* hyp,
2491                                      const TopoDS_Shape&             hypShape,
2492                                      set<TGeomID>&                   ignoreFaceIds)
2493 {
2494   TopExp_Explorer exp;
2495
2496   vector<TGeomID> ids = hyp->GetBndShapes();
2497   if ( hyp->IsToIgnoreShapes() ) // FACEs to ignore are given
2498   {
2499     for ( size_t ii = 0; ii < ids.size(); ++ii )
2500     {
2501       const TopoDS_Shape& s = getMeshDS()->IndexToShape( ids[ii] );
2502       if ( !s.IsNull() && s.ShapeType() == TopAbs_FACE )
2503         ignoreFaceIds.insert( ids[ii] );
2504     }
2505   }
2506   else // FACEs with layers are given
2507   {
2508     exp.Init( solid, TopAbs_FACE );
2509     for ( ; exp.More(); exp.Next() )
2510     {
2511       TGeomID faceInd = getMeshDS()->ShapeToIndex( exp.Current() );
2512       if ( find( ids.begin(), ids.end(), faceInd ) == ids.end() )
2513         ignoreFaceIds.insert( faceInd );
2514     }
2515   }
2516
2517   // ignore internal FACEs if inlets and outlets are specified
2518   if ( hyp->IsToIgnoreShapes() )
2519   {
2520     TopTools_IndexedDataMapOfShapeListOfShape solidsOfFace;
2521     TopExp::MapShapesAndAncestors( hypShape,
2522                                    TopAbs_FACE, TopAbs_SOLID, solidsOfFace);
2523
2524     for ( exp.Init( solid, TopAbs_FACE ); exp.More(); exp.Next() )
2525     {
2526       const TopoDS_Face& face = TopoDS::Face( exp.Current() );
2527       if ( SMESH_MesherHelper::NbAncestors( face, *_mesh, TopAbs_SOLID ) < 2 )
2528         continue;
2529
2530       int nbSolids = solidsOfFace.FindFromKey( face ).Extent();
2531       if ( nbSolids > 1 )
2532         ignoreFaceIds.insert( getMeshDS()->ShapeToIndex( face ));
2533     }
2534   }
2535 }
2536
2537 //================================================================================
2538 /*!
2539  * \brief Create the inner surface of the viscous layer and prepare data for infation
2540  */
2541 //================================================================================
2542
2543 bool _ViscousBuilder::makeLayer(_SolidData& data)
2544 {
2545   // get all sub-shapes to make layers on
2546   set<TGeomID> subIds, faceIds;
2547   subIds = data._noShrinkShapes;
2548   TopExp_Explorer exp( data._solid, TopAbs_FACE );
2549   for ( ; exp.More(); exp.Next() )
2550   {
2551     SMESH_subMesh* fSubM = _mesh->GetSubMesh( exp.Current() );
2552     if ( ! data._ignoreFaceIds.count( fSubM->GetId() ))
2553       faceIds.insert( fSubM->GetId() );
2554   }
2555
2556   // make a map to find new nodes on sub-shapes shared with other SOLID
2557   map< TGeomID, TNode2Edge* >::iterator s2ne;
2558   map< TGeomID, TopoDS_Shape >::iterator s2s = data._shrinkShape2Shape.begin();
2559   for (; s2s != data._shrinkShape2Shape.end(); ++s2s )
2560   {
2561     TGeomID shapeInd = s2s->first;
2562     for ( size_t i = 0; i < _sdVec.size(); ++i )
2563     {
2564       if ( _sdVec[i]._index == data._index ) continue;
2565       map< TGeomID, TopoDS_Shape >::iterator s2s2 = _sdVec[i]._shrinkShape2Shape.find( shapeInd );
2566       if ( s2s2 != _sdVec[i]._shrinkShape2Shape.end() &&
2567            *s2s == *s2s2 && !_sdVec[i]._n2eMap.empty() )
2568       {
2569         data._s2neMap.insert( make_pair( shapeInd, &_sdVec[i]._n2eMap ));
2570         break;
2571       }
2572     }
2573   }
2574
2575   // Create temporary faces and _LayerEdge's
2576
2577   dumpFunction(SMESH_Comment("makeLayers_")<<data._index);
2578
2579   data._stepSize = Precision::Infinite();
2580   data._stepSizeNodes[0] = 0;
2581
2582   SMESH_MesherHelper helper( *_mesh );
2583   helper.SetSubShape( data._solid );
2584   helper.SetElementsOnShape( true );
2585
2586   vector< const SMDS_MeshNode*> newNodes; // of a mesh face
2587   TNode2Edge::iterator n2e2;
2588
2589   // collect _LayerEdge's of shapes they are based on
2590   vector< _EdgesOnShape >& edgesByGeom = data._edgesOnShape;
2591   const int nbShapes = getMeshDS()->MaxShapeIndex();
2592   edgesByGeom.resize( nbShapes+1 );
2593
2594   // set data of _EdgesOnShape's
2595   if ( SMESH_subMesh* sm = _mesh->GetSubMesh( data._solid ))
2596   {
2597     SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(/*includeSelf=*/false);
2598     while ( smIt->more() )
2599     {
2600       sm = smIt->next();
2601       if ( sm->GetSubShape().ShapeType() == TopAbs_FACE &&
2602            !faceIds.count( sm->GetId() ))
2603         continue;
2604       setShapeData( edgesByGeom[ sm->GetId() ], sm, data );
2605     }
2606   }
2607   // make _LayerEdge's
2608   for ( set<TGeomID>::iterator id = faceIds.begin(); id != faceIds.end(); ++id )
2609   {
2610     const TopoDS_Face& F = TopoDS::Face( getMeshDS()->IndexToShape( *id ));
2611     SMESH_subMesh* sm = _mesh->GetSubMesh( F );
2612     SMESH_ProxyMesh::SubMesh* proxySub =
2613       data._proxyMesh->getFaceSubM( F, /*create=*/true);
2614
2615     SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
2616     if ( !smDS ) return error(SMESH_Comment("Not meshed face ") << *id, data._index );
2617
2618     SMDS_ElemIteratorPtr eIt = smDS->GetElements();
2619     while ( eIt->more() )
2620     {
2621       const SMDS_MeshElement* face = eIt->next();
2622       double          faceMaxCosin = -1;
2623       _LayerEdge*     maxCosinEdge = 0;
2624       int             nbDegenNodes = 0;
2625
2626       newNodes.resize( face->NbCornerNodes() );
2627       for ( size_t i = 0 ; i < newNodes.size(); ++i )
2628       {
2629         const SMDS_MeshNode* n = face->GetNode( i );
2630         const int      shapeID = n->getshapeId();
2631         const bool onDegenShap = helper.IsDegenShape( shapeID );
2632         const bool onDegenEdge = ( onDegenShap && n->GetPosition()->GetDim() == 1 );
2633         if ( onDegenShap )
2634         {
2635           if ( onDegenEdge )
2636           {
2637             // substitute n on a degenerated EDGE with a node on a corresponding VERTEX
2638             const TopoDS_Shape& E = getMeshDS()->IndexToShape( shapeID );
2639             TopoDS_Vertex       V = helper.IthVertex( 0, TopoDS::Edge( E ));
2640             if ( const SMDS_MeshNode* vN = SMESH_Algo::VertexNode( V, getMeshDS() )) {
2641               n = vN;
2642               nbDegenNodes++;
2643             }
2644           }
2645           else
2646           {
2647             nbDegenNodes++;
2648           }
2649         }
2650         TNode2Edge::iterator n2e = data._n2eMap.insert( make_pair( n, (_LayerEdge*)0 )).first;
2651         if ( !(*n2e).second )
2652         {
2653           // add a _LayerEdge
2654           _LayerEdge* edge = new _LayerEdge();
2655           edge->_nodes.push_back( n );
2656           n2e->second = edge;
2657           edgesByGeom[ shapeID ]._edges.push_back( edge );
2658           const bool noShrink = data._noShrinkShapes.count( shapeID );
2659
2660           SMESH_TNodeXYZ xyz( n );
2661
2662           // set edge data or find already refined _LayerEdge and get data from it
2663           if (( !noShrink                                                     ) &&
2664               ( n->GetPosition()->GetTypeOfPosition() != SMDS_TOP_FACE        ) &&
2665               (( s2ne = data._s2neMap.find( shapeID )) != data._s2neMap.end() ) &&
2666               (( n2e2 = (*s2ne).second->find( n )) != s2ne->second->end()     ))
2667           {
2668             _LayerEdge* foundEdge = (*n2e2).second;
2669             gp_XYZ        lastPos = edge->Copy( *foundEdge, edgesByGeom[ shapeID ], helper );
2670             foundEdge->_pos.push_back( lastPos );
2671             // location of the last node is modified and we restore it by foundEdge->_pos.back()
2672             const_cast< SMDS_MeshNode* >
2673               ( edge->_nodes.back() )->setXYZ( xyz.X(), xyz.Y(), xyz.Z() );
2674           }
2675           else
2676           {
2677             if ( !noShrink )
2678             {
2679               edge->_nodes.push_back( helper.AddNode( xyz.X(), xyz.Y(), xyz.Z() ));
2680             }
2681             if ( !setEdgeData( *edge, edgesByGeom[ shapeID ], helper, data ))
2682               return false;
2683
2684             if ( edge->_nodes.size() < 2 )
2685               edge->Block( data );
2686               //data._noShrinkShapes.insert( shapeID );
2687           }
2688           dumpMove(edge->_nodes.back());
2689
2690           if ( edge->_cosin > faceMaxCosin )
2691           {
2692             faceMaxCosin = edge->_cosin;
2693             maxCosinEdge = edge;
2694           }
2695         }
2696         newNodes[ i ] = n2e->second->_nodes.back();
2697
2698         if ( onDegenEdge )
2699           data._n2eMap.insert( make_pair( face->GetNode( i ), n2e->second ));
2700       }
2701       if ( newNodes.size() - nbDegenNodes < 2 )
2702         continue;
2703
2704       // create a temporary face
2705       const SMDS_MeshElement* newFace =
2706         new _TmpMeshFace( newNodes, --_tmpFaceID, face->GetShapeID(), face );
2707       proxySub->AddElement( newFace );
2708
2709       // compute inflation step size by min size of element on a convex surface
2710       if ( faceMaxCosin > theMinSmoothCosin )
2711         limitStepSize( data, face, maxCosinEdge );
2712
2713     } // loop on 2D elements on a FACE
2714   } // loop on FACEs of a SOLID to create _LayerEdge's
2715
2716
2717   // Set _LayerEdge::_neibors
2718   TNode2Edge::iterator n2e;
2719   for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
2720   {
2721     _EdgesOnShape& eos = data._edgesOnShape[iS];
2722     for ( size_t i = 0; i < eos._edges.size(); ++i )
2723     {
2724       _LayerEdge* edge = eos._edges[i];
2725       TIDSortedNodeSet nearNodes;
2726       SMDS_ElemIteratorPtr fIt = edge->_nodes[0]->GetInverseElementIterator(SMDSAbs_Face);
2727       while ( fIt->more() )
2728       {
2729         const SMDS_MeshElement* f = fIt->next();
2730         if ( !data._ignoreFaceIds.count( f->getshapeId() ))
2731           nearNodes.insert( f->begin_nodes(), f->end_nodes() );
2732       }
2733       nearNodes.erase( edge->_nodes[0] );
2734       edge->_neibors.reserve( nearNodes.size() );
2735       TIDSortedNodeSet::iterator node = nearNodes.begin();
2736       for ( ; node != nearNodes.end(); ++node )
2737         if (( n2e = data._n2eMap.find( *node )) != data._n2eMap.end() )
2738           edge->_neibors.push_back( n2e->second );
2739     }
2740   }
2741
2742   data._epsilon = 1e-7;
2743   if ( data._stepSize < 1. )
2744     data._epsilon *= data._stepSize;
2745
2746   if ( !findShapesToSmooth( data )) // _LayerEdge::_maxLen is computed here
2747     return false;
2748
2749   // limit data._stepSize depending on surface curvature and fill data._convexFaces
2750   limitStepSizeByCurvature( data ); // !!! it must be before node substitution in _Simplex
2751
2752   // Set target nodes into _Simplex and _LayerEdge's to _2NearEdges
2753   const SMDS_MeshNode* nn[2];
2754   for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
2755   {
2756     _EdgesOnShape& eos = data._edgesOnShape[iS];
2757     for ( size_t i = 0; i < eos._edges.size(); ++i )
2758     {
2759       _LayerEdge* edge = eos._edges[i];
2760       if ( edge->IsOnEdge() )
2761       {
2762         // get neighbor nodes
2763         bool hasData = ( edge->_2neibors->_edges[0] );
2764         if ( hasData ) // _LayerEdge is a copy of another one
2765         {
2766           nn[0] = edge->_2neibors->srcNode(0);
2767           nn[1] = edge->_2neibors->srcNode(1);
2768         }
2769         else if ( !findNeiborsOnEdge( edge, nn[0],nn[1], eos, data ))
2770         {
2771           return false;
2772         }
2773         // set neighbor _LayerEdge's
2774         for ( int j = 0; j < 2; ++j )
2775         {
2776           if (( n2e = data._n2eMap.find( nn[j] )) == data._n2eMap.end() )
2777             return error("_LayerEdge not found by src node", data._index);
2778           edge->_2neibors->_edges[j] = n2e->second;
2779         }
2780         if ( !hasData )
2781           edge->SetDataByNeighbors( nn[0], nn[1], eos, helper );
2782       }
2783
2784       for ( size_t j = 0; j < edge->_simplices.size(); ++j )
2785       {
2786         _Simplex& s = edge->_simplices[j];
2787         s._nNext = data._n2eMap[ s._nNext ]->_nodes.back();
2788         s._nPrev = data._n2eMap[ s._nPrev ]->_nodes.back();
2789       }
2790
2791       // For an _LayerEdge on a degenerated EDGE, copy some data from
2792       // a corresponding _LayerEdge on a VERTEX
2793       // (issue 52453, pb on a downloaded SampleCase2-Tet-netgen-mephisto.hdf)
2794       if ( helper.IsDegenShape( edge->_nodes[0]->getshapeId() ))
2795       {
2796         // Generally we should not get here
2797         if ( eos.ShapeType() != TopAbs_EDGE )
2798           continue;
2799         TopoDS_Vertex V = helper.IthVertex( 0, TopoDS::Edge( eos._shape ));
2800         const SMDS_MeshNode* vN = SMESH_Algo::VertexNode( V, getMeshDS() );
2801         if (( n2e = data._n2eMap.find( vN )) == data._n2eMap.end() )
2802           continue;
2803         const _LayerEdge* vEdge = n2e->second;
2804         edge->_normal    = vEdge->_normal;
2805         edge->_lenFactor = vEdge->_lenFactor;
2806         edge->_cosin     = vEdge->_cosin;
2807       }
2808
2809     } // loop on data._edgesOnShape._edges
2810   } // loop on data._edgesOnShape
2811
2812   // fix _LayerEdge::_2neibors on EDGEs to smooth
2813   // map< TGeomID,Handle(Geom_Curve)>::iterator e2c = data._edge2curve.begin();
2814   // for ( ; e2c != data._edge2curve.end(); ++e2c )
2815   //   if ( !e2c->second.IsNull() )
2816   //   {
2817   //     if ( _EdgesOnShape* eos = data.GetShapeEdges( e2c->first ))
2818   //       data.Sort2NeiborsOnEdge( eos->_edges );
2819   //   }
2820
2821   dumpFunctionEnd();
2822   return true;
2823 }
2824
2825 //================================================================================
2826 /*!
2827  * \brief Compute inflation step size by min size of element on a convex surface
2828  */
2829 //================================================================================
2830
2831 void _ViscousBuilder::limitStepSize( _SolidData&             data,
2832                                      const SMDS_MeshElement* face,
2833                                      const _LayerEdge*       maxCosinEdge )
2834 {
2835   int iN = 0;
2836   double minSize = 10 * data._stepSize;
2837   const int nbNodes = face->NbCornerNodes();
2838   for ( int i = 0; i < nbNodes; ++i )
2839   {
2840     const SMDS_MeshNode* nextN = face->GetNode( SMESH_MesherHelper::WrapIndex( i+1, nbNodes ));
2841     const SMDS_MeshNode*  curN = face->GetNode( i );
2842     if ( nextN->GetPosition()->GetTypeOfPosition() == SMDS_TOP_FACE ||
2843          curN-> GetPosition()->GetTypeOfPosition() == SMDS_TOP_FACE )
2844     {
2845       double dist = SMESH_TNodeXYZ( curN ).Distance( nextN );
2846       if ( dist < minSize )
2847         minSize = dist, iN = i;
2848     }
2849   }
2850   double newStep = 0.8 * minSize / maxCosinEdge->_lenFactor;
2851   if ( newStep < data._stepSize )
2852   {
2853     data._stepSize = newStep;
2854     data._stepSizeCoeff = 0.8 / maxCosinEdge->_lenFactor;
2855     data._stepSizeNodes[0] = face->GetNode( iN );
2856     data._stepSizeNodes[1] = face->GetNode( SMESH_MesherHelper::WrapIndex( iN+1, nbNodes ));
2857   }
2858 }
2859
2860 //================================================================================
2861 /*!
2862  * \brief Compute inflation step size by min size of element on a convex surface
2863  */
2864 //================================================================================
2865
2866 void _ViscousBuilder::limitStepSize( _SolidData& data, const double minSize )
2867 {
2868   if ( minSize < data._stepSize )
2869   {
2870     data._stepSize = minSize;
2871     if ( data._stepSizeNodes[0] )
2872     {
2873       double dist =
2874         SMESH_TNodeXYZ(data._stepSizeNodes[0]).Distance(data._stepSizeNodes[1]);
2875       data._stepSizeCoeff = data._stepSize / dist;
2876     }
2877   }
2878 }
2879
2880 //================================================================================
2881 /*!
2882  * \brief Limit data._stepSize by evaluating curvature of shapes and fill data._convexFaces
2883  */
2884 //================================================================================
2885
2886 void _ViscousBuilder::limitStepSizeByCurvature( _SolidData& data )
2887 {
2888   SMESH_MesherHelper helper( *_mesh );
2889
2890   BRepLProp_SLProps surfProp( 2, 1e-6 );
2891   data._convexFaces.clear();
2892
2893   for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
2894   {
2895     _EdgesOnShape& eof = data._edgesOnShape[iS];
2896     if ( eof.ShapeType() != TopAbs_FACE ||
2897          data._ignoreFaceIds.count( eof._shapeID ))
2898       continue;
2899
2900     TopoDS_Face        F = TopoDS::Face( eof._shape );
2901     const TGeomID faceID = eof._shapeID;
2902
2903     BRepAdaptor_Surface surface( F, false );
2904     surfProp.SetSurface( surface );
2905
2906     _ConvexFace cnvFace;
2907     cnvFace._face = F;
2908     cnvFace._normalsFixed = false;
2909     cnvFace._isTooCurved = false;
2910
2911     double maxCurvature = cnvFace.GetMaxCurvature( data, eof, surfProp, helper );
2912     if ( maxCurvature > 0 )
2913     {
2914       limitStepSize( data, 0.9 / maxCurvature );
2915       findEdgesToUpdateNormalNearConvexFace( cnvFace, data, helper );
2916     }
2917     if ( !cnvFace._isTooCurved ) continue;
2918
2919     _ConvexFace & convFace =
2920       data._convexFaces.insert( make_pair( faceID, cnvFace )).first->second;
2921
2922     // skip a closed surface (data._convexFaces is useful anyway)
2923     bool isClosedF = false;
2924     helper.SetSubShape( F );
2925     if ( helper.HasRealSeam() )
2926     {
2927       // in the closed surface there must be a closed EDGE
2928       for ( TopExp_Explorer eIt( F, TopAbs_EDGE ); eIt.More() && !isClosedF; eIt.Next() )
2929         isClosedF = helper.IsClosedEdge( TopoDS::Edge( eIt.Current() ));
2930     }
2931     if ( isClosedF )
2932     {
2933       // limit _LayerEdge::_maxLen on the FACE
2934       const double oriFactor    = ( F.Orientation() == TopAbs_REVERSED ? +1. : -1. );
2935       const double minCurvature =
2936         1. / ( eof._hyp.GetTotalThickness() * ( 1 + theThickToIntersection ));
2937       map< TGeomID, _EdgesOnShape* >::iterator id2eos = cnvFace._subIdToEOS.find( faceID );
2938       if ( id2eos != cnvFace._subIdToEOS.end() )
2939       {
2940         _EdgesOnShape& eos = * id2eos->second;
2941         for ( size_t i = 0; i < eos._edges.size(); ++i )
2942         {
2943           _LayerEdge* ledge = eos._edges[ i ];
2944           gp_XY uv = helper.GetNodeUV( F, ledge->_nodes[0] );
2945           surfProp.SetParameters( uv.X(), uv.Y() );
2946           if ( surfProp.IsCurvatureDefined() )
2947           {
2948             double curvature = Max( surfProp.MaxCurvature() * oriFactor,
2949                                     surfProp.MinCurvature() * oriFactor );
2950             if ( curvature > minCurvature )
2951               ledge->SetMaxLen( Min( ledge->_maxLen, 1. / curvature ));
2952           }
2953         }
2954       }
2955       continue;
2956     }
2957
2958     // Fill _ConvexFace::_simplexTestEdges. These _LayerEdge's are used to detect
2959     // prism distortion.
2960     map< TGeomID, _EdgesOnShape* >::iterator id2eos = convFace._subIdToEOS.find( faceID );
2961     if ( id2eos != convFace._subIdToEOS.end() && !id2eos->second->_edges.empty() )
2962     {
2963       // there are _LayerEdge's on the FACE it-self;
2964       // select _LayerEdge's near EDGEs
2965       _EdgesOnShape& eos = * id2eos->second;
2966       for ( size_t i = 0; i < eos._edges.size(); ++i )
2967       {
2968         _LayerEdge* ledge = eos._edges[ i ];
2969         for ( size_t j = 0; j < ledge->_simplices.size(); ++j )
2970           if ( ledge->_simplices[j]._nNext->GetPosition()->GetDim() < 2 )
2971           {
2972             // do not select _LayerEdge's neighboring sharp EDGEs
2973             bool sharpNbr = false;
2974             for ( size_t iN = 0; iN < ledge->_neibors.size()  && !sharpNbr; ++iN )
2975               sharpNbr = ( ledge->_neibors[iN]->_cosin > theMinSmoothCosin );
2976             if ( !sharpNbr )
2977               convFace._simplexTestEdges.push_back( ledge );
2978             break;
2979           }
2980       }
2981     }
2982     else
2983     {
2984       // where there are no _LayerEdge's on a _ConvexFace,
2985       // as e.g. on a fillet surface with no internal nodes - issue 22580,
2986       // so that collision of viscous internal faces is not detected by check of
2987       // intersection of _LayerEdge's with the viscous internal faces.
2988
2989       set< const SMDS_MeshNode* > usedNodes;
2990
2991       // look for _LayerEdge's with null _sWOL
2992       id2eos = convFace._subIdToEOS.begin();
2993       for ( ; id2eos != convFace._subIdToEOS.end(); ++id2eos )
2994       {
2995         _EdgesOnShape& eos = * id2eos->second;
2996         if ( !eos._sWOL.IsNull() )
2997           continue;
2998         for ( size_t i = 0; i < eos._edges.size(); ++i )
2999         {
3000           _LayerEdge* ledge = eos._edges[ i ];
3001           const SMDS_MeshNode* srcNode = ledge->_nodes[0];
3002           if ( !usedNodes.insert( srcNode ).second ) continue;
3003
3004           for ( size_t i = 0; i < ledge->_simplices.size(); ++i )
3005           {
3006             usedNodes.insert( ledge->_simplices[i]._nPrev );
3007             usedNodes.insert( ledge->_simplices[i]._nNext );
3008           }
3009           convFace._simplexTestEdges.push_back( ledge );
3010         }
3011       }
3012     }
3013   } // loop on FACEs of data._solid
3014 }
3015
3016 //================================================================================
3017 /*!
3018  * \brief Detect shapes (and _LayerEdge's on them) to smooth
3019  */
3020 //================================================================================
3021
3022 bool _ViscousBuilder::findShapesToSmooth( _SolidData& data )
3023 {
3024   // define allowed thickness
3025   computeGeomSize( data ); // compute data._geomSize and _LayerEdge::_maxLen
3026
3027
3028   // Find shapes needing smoothing; such a shape has _LayerEdge._normal on it's
3029   // boundary inclined to the shape at a sharp angle
3030
3031   TopTools_MapOfShape edgesOfSmooFaces;
3032   SMESH_MesherHelper helper( *_mesh );
3033   bool ok = true;
3034
3035   vector< _EdgesOnShape >& edgesByGeom = data._edgesOnShape;
3036   data._nbShapesToSmooth = 0;
3037
3038   for ( size_t iS = 0; iS < edgesByGeom.size(); ++iS ) // check FACEs
3039   {
3040     _EdgesOnShape& eos = edgesByGeom[iS];
3041     eos._toSmooth = false;
3042     if ( eos._edges.empty() || eos.ShapeType() != TopAbs_FACE )
3043       continue;
3044
3045     double tgtThick = eos._hyp.GetTotalThickness();
3046     SMESH_subMeshIteratorPtr subIt = eos._subMesh->getDependsOnIterator(/*includeSelf=*/false );
3047     while ( subIt->more() && !eos._toSmooth )
3048     {
3049       TGeomID iSub = subIt->next()->GetId();
3050       const vector<_LayerEdge*>& eSub = edgesByGeom[ iSub ]._edges;
3051       if ( eSub.empty() ) continue;
3052
3053       double faceSize;
3054       for ( size_t i = 0; i < eSub.size() && !eos._toSmooth; ++i )
3055         if ( eSub[i]->_cosin > theMinSmoothCosin )
3056         {
3057           SMDS_ElemIteratorPtr fIt = eSub[i]->_nodes[0]->GetInverseElementIterator(SMDSAbs_Face);
3058           while ( fIt->more() && !eos._toSmooth )
3059           {
3060             const SMDS_MeshElement* face = fIt->next();
3061             if ( face->getshapeId() == eos._shapeID &&
3062                  getDistFromEdge( face, eSub[i]->_nodes[0], faceSize ))
3063             {
3064               eos._toSmooth = needSmoothing( eSub[i]->_cosin,
3065                                              tgtThick * eSub[i]->_lenFactor,
3066                                              faceSize);
3067             }
3068           }
3069         }
3070     }
3071     if ( eos._toSmooth )
3072     {
3073       for ( TopExp_Explorer eExp( edgesByGeom[iS]._shape, TopAbs_EDGE ); eExp.More(); eExp.Next() )
3074         edgesOfSmooFaces.Add( eExp.Current() );
3075
3076       data.PrepareEdgesToSmoothOnFace( &edgesByGeom[iS], /*substituteSrcNodes=*/false );
3077     }
3078     data._nbShapesToSmooth += eos._toSmooth;
3079
3080   }  // check FACEs
3081
3082   for ( size_t iS = 0; iS < edgesByGeom.size(); ++iS ) // check EDGEs
3083   {
3084     _EdgesOnShape& eos = edgesByGeom[iS];
3085     eos._edgeSmoother = NULL;
3086     if ( eos._edges.empty() || eos.ShapeType() != TopAbs_EDGE ) continue;
3087     if ( !eos._hyp.ToSmooth() ) continue;
3088
3089     const TopoDS_Edge& E = TopoDS::Edge( edgesByGeom[iS]._shape );
3090     if ( SMESH_Algo::isDegenerated( E ) || !edgesOfSmooFaces.Contains( E ))
3091       continue;
3092
3093     double tgtThick = eos._hyp.GetTotalThickness();
3094     for ( TopoDS_Iterator vIt( E ); vIt.More() && !eos._toSmooth; vIt.Next() )
3095     {
3096       TGeomID iV = getMeshDS()->ShapeToIndex( vIt.Value() );
3097       vector<_LayerEdge*>& eV = edgesByGeom[ iV ]._edges;
3098       if ( eV.empty() || eV[0]->Is( _LayerEdge::MULTI_NORMAL )) continue;
3099       gp_Vec  eDir    = getEdgeDir( E, TopoDS::Vertex( vIt.Value() ));
3100       double angle    = eDir.Angle( eV[0]->_normal );
3101       double cosin    = Cos( angle );
3102       double cosinAbs = Abs( cosin );
3103       if ( cosinAbs > theMinSmoothCosin )
3104       {
3105         // always smooth analytic EDGEs
3106         Handle(Geom_Curve) curve = _Smoother1D::CurveForSmooth( E, eos, helper );
3107         eos._toSmooth = ! curve.IsNull();
3108
3109         // compare tgtThick with the length of an end segment
3110         SMDS_ElemIteratorPtr eIt = eV[0]->_nodes[0]->GetInverseElementIterator(SMDSAbs_Edge);
3111         while ( eIt->more() && !eos._toSmooth )
3112         {
3113           const SMDS_MeshElement* endSeg = eIt->next();
3114           if ( endSeg->getshapeId() == (int) iS )
3115           {
3116             double segLen =
3117               SMESH_TNodeXYZ( endSeg->GetNode( 0 )).Distance( endSeg->GetNode( 1 ));
3118             eos._toSmooth = needSmoothing( cosinAbs, tgtThick * eV[0]->_lenFactor, segLen );
3119           }
3120         }
3121         if ( eos._toSmooth )
3122         {
3123           eos._edgeSmoother = new _Smoother1D( curve, eos );
3124
3125           // for ( size_t i = 0; i < eos._edges.size(); ++i )
3126           //   eos._edges[i]->Set( _LayerEdge::TO_SMOOTH );
3127         }
3128       }
3129     }
3130     data._nbShapesToSmooth += eos._toSmooth;
3131
3132   } // check EDGEs
3133
3134   // Reset _cosin if no smooth is allowed by the user
3135   for ( size_t iS = 0; iS < edgesByGeom.size(); ++iS )
3136   {
3137     _EdgesOnShape& eos = edgesByGeom[iS];
3138     if ( eos._edges.empty() ) continue;
3139
3140     if ( !eos._hyp.ToSmooth() )
3141       for ( size_t i = 0; i < eos._edges.size(); ++i )
3142         //eos._edges[i]->SetCosin( 0 ); // keep _cosin to use in limitMaxLenByCurvature()
3143         eos._edges[i]->_lenFactor = 1;
3144   }
3145
3146
3147   // Fill _eosC1 to make that C1 FACEs and EDGEs between them to be smoothed as a whole
3148
3149   TopTools_MapOfShape c1VV;
3150
3151   for ( size_t iS = 0; iS < edgesByGeom.size(); ++iS ) // check FACEs
3152   {
3153     _EdgesOnShape& eos = edgesByGeom[iS];
3154     if ( eos._edges.empty() ||
3155          eos.ShapeType() != TopAbs_FACE ||
3156          !eos._toSmooth )
3157       continue;
3158
3159     // check EDGEs of a FACE
3160     TopTools_MapOfShape checkedEE, allVV;
3161     list< SMESH_subMesh* > smQueue( 1, eos._subMesh ); // sm of FACEs
3162     while ( !smQueue.empty() )
3163     {
3164       SMESH_subMesh* sm = smQueue.front();
3165       smQueue.pop_front();
3166       SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(/*includeSelf=*/false);
3167       while ( smIt->more() )
3168       {
3169         sm = smIt->next();
3170         if ( sm->GetSubShape().ShapeType() == TopAbs_VERTEX )
3171           allVV.Add( sm->GetSubShape() );
3172         if ( sm->GetSubShape().ShapeType() != TopAbs_EDGE ||
3173              !checkedEE.Add( sm->GetSubShape() ))
3174           continue;
3175
3176         _EdgesOnShape*      eoe = data.GetShapeEdges( sm->GetId() );
3177         vector<_LayerEdge*>& eE = eoe->_edges;
3178         if ( eE.empty() || !eoe->_sWOL.IsNull() )
3179           continue;
3180
3181         bool isC1 = true; // check continuity along an EDGE
3182         for ( size_t i = 0; i < eE.size() && isC1; ++i )
3183           isC1 = ( Abs( eE[i]->_cosin ) < theMinSmoothCosin );
3184         if ( !isC1 )
3185           continue;
3186
3187         // check that mesh faces are C1 as well
3188         {
3189           gp_XYZ norm1, norm2;
3190           const SMDS_MeshNode*   n = eE[ eE.size() / 2 ]->_nodes[0];
3191           SMDS_ElemIteratorPtr fIt = n->GetInverseElementIterator(SMDSAbs_Face);
3192           if ( !SMESH_MeshAlgos::FaceNormal( fIt->next(), norm1, /*normalized=*/true ))
3193             continue;
3194           while ( fIt->more() && isC1 )
3195             isC1 = ( SMESH_MeshAlgos::FaceNormal( fIt->next(), norm2, /*normalized=*/true ) &&
3196                      Abs( norm1 * norm2 ) >= ( 1. - theMinSmoothCosin ));
3197           if ( !isC1 )
3198             continue;
3199         }
3200
3201         // add the EDGE and an adjacent FACE to _eosC1
3202         PShapeIteratorPtr fIt = helper.GetAncestors( sm->GetSubShape(), *_mesh, TopAbs_FACE );
3203         while ( const TopoDS_Shape* face = fIt->next() )
3204         {
3205           _EdgesOnShape* eof = data.GetShapeEdges( *face );
3206           if ( !eof ) continue; // other solid
3207           if ( eos._shapeID == eof->_shapeID ) continue;
3208           if ( !eos.HasC1( eof ))
3209           {
3210             // check the FACEs
3211             eos._eosC1.push_back( eof );
3212             eof->_toSmooth = false;
3213             data.PrepareEdgesToSmoothOnFace( eof, /*substituteSrcNodes=*/false );
3214             smQueue.push_back( eof->_subMesh );
3215           }
3216           if ( !eos.HasC1( eoe ))
3217           {
3218             eos._eosC1.push_back( eoe );
3219             eoe->_toSmooth = false;
3220             data.PrepareEdgesToSmoothOnFace( eoe, /*substituteSrcNodes=*/false );
3221           }
3222         }
3223       }
3224     }
3225     if ( eos._eosC1.empty() )
3226       continue;
3227
3228     // check VERTEXes of C1 FACEs
3229     TopTools_MapIteratorOfMapOfShape vIt( allVV );
3230     for ( ; vIt.More(); vIt.Next() )
3231     {
3232       _EdgesOnShape* eov = data.GetShapeEdges( vIt.Key() );
3233       if ( !eov || eov->_edges.empty() || !eov->_sWOL.IsNull() )
3234         continue;
3235
3236       bool isC1 = true; // check if all adjacent FACEs are in eos._eosC1
3237       PShapeIteratorPtr fIt = helper.GetAncestors( vIt.Key(), *_mesh, TopAbs_FACE );
3238       while ( const TopoDS_Shape* face = fIt->next() )
3239       {
3240         _EdgesOnShape* eof = data.GetShapeEdges( *face );
3241         if ( !eof ) continue; // other solid
3242         isC1 = ( face->IsSame( eos._shape ) || eos.HasC1( eof ));
3243         if ( !isC1 )
3244           break;
3245       }
3246       if ( isC1 )
3247       {
3248         eos._eosC1.push_back( eov );
3249         data.PrepareEdgesToSmoothOnFace( eov, /*substituteSrcNodes=*/false );
3250         c1VV.Add( eov->_shape );
3251       }
3252     }
3253
3254   } // fill _eosC1 of FACEs
3255
3256
3257   // Find C1 EDGEs
3258
3259   vector< pair< _EdgesOnShape*, gp_XYZ > > dirOfEdges;
3260
3261   for ( size_t iS = 0; iS < edgesByGeom.size(); ++iS ) // check VERTEXes
3262   {
3263     _EdgesOnShape& eov = edgesByGeom[iS];
3264     if ( eov._edges.empty() ||
3265          eov.ShapeType() != TopAbs_VERTEX ||
3266          c1VV.Contains( eov._shape ))
3267       continue;
3268     const TopoDS_Vertex& V = TopoDS::Vertex( eov._shape );
3269
3270     // get directions of surrounding EDGEs
3271     dirOfEdges.clear();
3272     PShapeIteratorPtr fIt = helper.GetAncestors( eov._shape, *_mesh, TopAbs_EDGE );
3273     while ( const TopoDS_Shape* e = fIt->next() )
3274     {
3275       _EdgesOnShape* eoe = data.GetShapeEdges( *e );
3276       if ( !eoe ) continue; // other solid
3277       gp_XYZ eDir = getEdgeDir( TopoDS::Edge( *e ), V );
3278       if ( !Precision::IsInfinite( eDir.X() ))
3279         dirOfEdges.push_back( make_pair( eoe, eDir.Normalized() ));
3280     }
3281
3282     // find EDGEs with C1 directions
3283     for ( size_t i = 0; i < dirOfEdges.size(); ++i )
3284       for ( size_t j = i+1; j < dirOfEdges.size(); ++j )
3285         if ( dirOfEdges[i].first && dirOfEdges[j].first )
3286         {
3287           double dot = dirOfEdges[i].second * dirOfEdges[j].second;
3288           bool isC1 = ( dot < - ( 1. - theMinSmoothCosin ));
3289           if ( isC1 )
3290           {
3291             double maxEdgeLen = 3 * Min( eov._edges[0]->_maxLen, eov._hyp.GetTotalThickness() );
3292             double eLen1 = SMESH_Algo::EdgeLength( TopoDS::Edge( dirOfEdges[i].first->_shape ));
3293             double eLen2 = SMESH_Algo::EdgeLength( TopoDS::Edge( dirOfEdges[j].first->_shape ));
3294             if ( eLen1 < maxEdgeLen ) eov._eosC1.push_back( dirOfEdges[i].first );
3295             if ( eLen2 < maxEdgeLen ) eov._eosC1.push_back( dirOfEdges[j].first );
3296             dirOfEdges[i].first = 0;
3297             dirOfEdges[j].first = 0;
3298           }
3299         }
3300   } // fill _eosC1 of VERTEXes
3301
3302
3303
3304   return ok;
3305 }
3306
3307 //================================================================================
3308 /*!
3309  * \brief initialize data of _EdgesOnShape
3310  */
3311 //================================================================================
3312
3313 void _ViscousBuilder::setShapeData( _EdgesOnShape& eos,
3314                                     SMESH_subMesh* sm,
3315                                     _SolidData&    data )
3316 {
3317   if ( !eos._shape.IsNull() ||
3318        sm->GetSubShape().ShapeType() == TopAbs_WIRE )
3319     return;
3320
3321   SMESH_MesherHelper helper( *_mesh );
3322
3323   eos._subMesh = sm;
3324   eos._shapeID = sm->GetId();
3325   eos._shape   = sm->GetSubShape();
3326   if ( eos.ShapeType() == TopAbs_FACE )
3327     eos._shape.Orientation( helper.GetSubShapeOri( data._solid, eos._shape ));
3328   eos._toSmooth = false;
3329   eos._data = &data;
3330
3331   // set _SWOL
3332   map< TGeomID, TopoDS_Shape >::const_iterator s2s =
3333     data._shrinkShape2Shape.find( eos._shapeID );
3334   if ( s2s != data._shrinkShape2Shape.end() )
3335     eos._sWOL = s2s->second;
3336
3337   eos._isRegularSWOL = true;
3338   if ( eos.SWOLType() == TopAbs_FACE )
3339   {
3340     const TopoDS_Face& F = TopoDS::Face( eos._sWOL );
3341     Handle(ShapeAnalysis_Surface) surface = helper.GetSurface( F );
3342     eos._isRegularSWOL = ( ! surface->HasSingularities( 1e-7 ));
3343   }
3344
3345   // set _hyp
3346   if ( data._hyps.size() == 1 )
3347   {
3348     eos._hyp = data._hyps.back();
3349   }
3350   else
3351   {
3352     // compute average StdMeshers_ViscousLayers parameters
3353     map< TGeomID, const StdMeshers_ViscousLayers* >::iterator f2hyp;
3354     if ( eos.ShapeType() == TopAbs_FACE )
3355     {
3356       if (( f2hyp = data._face2hyp.find( eos._shapeID )) != data._face2hyp.end() )
3357         eos._hyp = f2hyp->second;
3358     }
3359     else
3360     {
3361       PShapeIteratorPtr fIt = helper.GetAncestors( eos._shape, *_mesh, TopAbs_FACE );
3362       while ( const TopoDS_Shape* face = fIt->next() )
3363       {
3364         TGeomID faceID = getMeshDS()->ShapeToIndex( *face );
3365         if (( f2hyp = data._face2hyp.find( faceID )) != data._face2hyp.end() )
3366           eos._hyp.Add( f2hyp->second );
3367       }
3368     }
3369   }
3370
3371   // set _faceNormals
3372   if ( ! eos._hyp.UseSurfaceNormal() )
3373   {
3374     if ( eos.ShapeType() == TopAbs_FACE ) // get normals to elements on a FACE
3375     {
3376       SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
3377       if ( !smDS ) return;
3378       eos._faceNormals.reserve( smDS->NbElements() );
3379
3380       double oriFactor = helper.IsReversedSubMesh( TopoDS::Face( eos._shape )) ? 1.: -1.;
3381       SMDS_ElemIteratorPtr eIt = smDS->GetElements();
3382       for ( ; eIt->more(); )
3383       {
3384         const SMDS_MeshElement* face = eIt->next();
3385         gp_XYZ&                 norm = eos._faceNormals[face];
3386         if ( !SMESH_MeshAlgos::FaceNormal( face, norm, /*normalized=*/true ))
3387           norm.SetCoord( 0,0,0 );
3388         norm *= oriFactor;
3389       }
3390     }
3391     else // find EOS of adjacent FACEs
3392     {
3393       PShapeIteratorPtr fIt = helper.GetAncestors( eos._shape, *_mesh, TopAbs_FACE );
3394       while ( const TopoDS_Shape* face = fIt->next() )
3395       {
3396         TGeomID faceID = getMeshDS()->ShapeToIndex( *face );
3397         eos._faceEOS.push_back( & data._edgesOnShape[ faceID ]);
3398         if ( eos._faceEOS.back()->_shape.IsNull() )
3399           // avoid using uninitialised _shapeID in GetNormal()
3400           eos._faceEOS.back()->_shapeID = faceID;
3401       }
3402     }
3403   }
3404 }
3405
3406 //================================================================================
3407 /*!
3408  * \brief Returns normal of a face
3409  */
3410 //================================================================================
3411
3412 bool _EdgesOnShape::GetNormal( const SMDS_MeshElement* face, gp_Vec& norm )
3413 {
3414   bool ok = false;
3415   _EdgesOnShape* eos = 0;
3416
3417   if ( face->getshapeId() == _shapeID )
3418   {
3419     eos = this;
3420   }
3421   else
3422   {
3423     for ( size_t iF = 0; iF < _faceEOS.size() && !eos; ++iF )
3424       if ( face->getshapeId() == _faceEOS[ iF ]->_shapeID )
3425         eos = _faceEOS[ iF ];
3426   }
3427
3428   if (( eos ) &&
3429       ( ok = ( eos->_faceNormals.count( face ) )))
3430   {
3431     norm = eos->_faceNormals[ face ];
3432   }
3433   else if ( !eos )
3434   {
3435     debugMsg( "_EdgesOnShape::Normal() failed for face "<<face->GetID()
3436               << " on _shape #" << _shapeID );
3437   }
3438   return ok;
3439 }
3440
3441
3442 //================================================================================
3443 /*!
3444  * \brief Set data of _LayerEdge needed for smoothing
3445  */
3446 //================================================================================
3447
3448 bool _ViscousBuilder::setEdgeData(_LayerEdge&         edge,
3449                                   _EdgesOnShape&      eos,
3450                                   SMESH_MesherHelper& helper,
3451                                   _SolidData&         data)
3452 {
3453   const SMDS_MeshNode* node = edge._nodes[0]; // source node
3454
3455   edge._len       = 0;
3456   edge._maxLen    = Precision::Infinite();
3457   edge._minAngle  = 0;
3458   edge._2neibors  = 0;
3459   edge._curvature = 0;
3460   edge._flags     = 0;
3461
3462   // --------------------------
3463   // Compute _normal and _cosin
3464   // --------------------------
3465
3466   edge._cosin     = 0;
3467   edge._lenFactor = 1.;
3468   edge._normal.SetCoord(0,0,0);
3469   _Simplex::GetSimplices( node, edge._simplices, data._ignoreFaceIds, &data );
3470
3471   int totalNbFaces = 0;
3472   TopoDS_Face F;
3473   std::pair< TopoDS_Face, gp_XYZ > face2Norm[20];
3474   gp_Vec geomNorm;
3475   bool normOK = true;
3476
3477   const bool onShrinkShape = !eos._sWOL.IsNull();
3478   const bool useGeometry   = (( eos._hyp.UseSurfaceNormal() ) ||
3479                               ( eos.ShapeType() != TopAbs_FACE /*&& !onShrinkShape*/ ));
3480
3481   // get geom FACEs the node lies on
3482   //if ( useGeometry )
3483   {
3484     set<TGeomID> faceIds;
3485     if  ( eos.ShapeType() == TopAbs_FACE )
3486     {
3487       faceIds.insert( eos._shapeID );
3488     }
3489     else
3490     {
3491       SMDS_ElemIteratorPtr fIt = node->GetInverseElementIterator(SMDSAbs_Face);
3492       while ( fIt->more() )
3493         faceIds.insert( fIt->next()->getshapeId() );
3494     }
3495     set<TGeomID>::iterator id = faceIds.begin();
3496     for ( ; id != faceIds.end(); ++id )
3497     {
3498       const TopoDS_Shape& s = getMeshDS()->IndexToShape( *id );
3499       if ( s.IsNull() || s.ShapeType() != TopAbs_FACE || data._ignoreFaceIds.count( *id ))
3500         continue;
3501       F = TopoDS::Face( s );
3502       face2Norm[ totalNbFaces ].first = F;
3503       totalNbFaces++;
3504     }
3505   }
3506
3507   // find _normal
3508   bool fromVonF = false;
3509   if ( useGeometry )
3510   {
3511     fromVonF = ( eos.ShapeType() == TopAbs_VERTEX &&
3512                  eos.SWOLType()  == TopAbs_FACE  &&
3513                  totalNbFaces > 1 );
3514
3515     if ( onShrinkShape && !fromVonF ) // one of faces the node is on has no layers
3516     {
3517       if ( eos.SWOLType() == TopAbs_EDGE )
3518       {
3519         // inflate from VERTEX along EDGE
3520         edge._normal = getEdgeDir( TopoDS::Edge( eos._sWOL ), TopoDS::Vertex( eos._shape ));
3521       }
3522       else if ( eos.ShapeType() == TopAbs_VERTEX )
3523       {
3524         // inflate from VERTEX along FACE
3525         edge._normal = getFaceDir( TopoDS::Face( eos._sWOL ), TopoDS::Vertex( eos._shape ),
3526                                    node, helper, normOK, &edge._cosin);
3527       }
3528       else
3529       {
3530         // inflate from EDGE along FACE
3531         edge._normal = getFaceDir( TopoDS::Face( eos._sWOL ), TopoDS::Edge( eos._shape ),
3532                                    node, helper, normOK);
3533       }
3534     }
3535     else // layers are on all FACEs of SOLID the node is on (or fromVonF)
3536     {
3537       if ( fromVonF )
3538         face2Norm[ totalNbFaces++ ].first = TopoDS::Face( eos._sWOL );
3539
3540       int nbOkNorms = 0;
3541       for ( int iF = totalNbFaces - 1; iF >= 0; --iF )
3542       {
3543         F = face2Norm[ iF ].first;
3544         geomNorm = getFaceNormal( node, F, helper, normOK );
3545         if ( !normOK ) continue;
3546         nbOkNorms++;
3547
3548         if ( helper.GetSubShapeOri( data._solid, F ) != TopAbs_REVERSED )
3549           geomNorm.Reverse();
3550         face2Norm[ iF ].second = geomNorm.XYZ();
3551         edge._normal += geomNorm.XYZ();
3552       }
3553       if ( nbOkNorms == 0 )
3554         return error(SMESH_Comment("Can't get normal to node ") << node->GetID(), data._index);
3555
3556       if ( totalNbFaces >= 3 )
3557       {
3558         edge._normal = getNormalByOffset( &edge, face2Norm, totalNbFaces, fromVonF );
3559       }
3560
3561       if ( edge._normal.Modulus() < 1e-3 && nbOkNorms > 1 )
3562       {
3563         // opposite normals, re-get normals at shifted positions (IPAL 52426)
3564         edge._normal.SetCoord( 0,0,0 );
3565         for ( int iF = 0; iF < totalNbFaces - fromVonF; ++iF )
3566         {
3567           const TopoDS_Face& F = face2Norm[iF].first;
3568           geomNorm = getFaceNormal( node, F, helper, normOK, /*shiftInside=*/true );
3569           if ( helper.GetSubShapeOri( data._solid, F ) != TopAbs_REVERSED )
3570             geomNorm.Reverse();
3571           if ( normOK )
3572             face2Norm[ iF ].second = geomNorm.XYZ();
3573           edge._normal += face2Norm[ iF ].second;
3574         }
3575       }
3576     }
3577   }
3578   else // !useGeometry - get _normal using surrounding mesh faces
3579   {
3580     edge._normal = getWeigthedNormal( &edge );
3581
3582     // set<TGeomID> faceIds;
3583     //
3584     // SMDS_ElemIteratorPtr fIt = node->GetInverseElementIterator(SMDSAbs_Face);
3585     // while ( fIt->more() )
3586     // {
3587     //   const SMDS_MeshElement* face = fIt->next();
3588     //   if ( eos.GetNormal( face, geomNorm ))
3589     //   {
3590     //     if ( onShrinkShape && !faceIds.insert( face->getshapeId() ).second )
3591     //       continue; // use only one mesh face on FACE
3592     //     edge._normal += geomNorm.XYZ();
3593     //     totalNbFaces++;
3594     //   }
3595     // }
3596   }
3597
3598   // compute _cosin
3599   //if ( eos._hyp.UseSurfaceNormal() )
3600   {
3601     switch ( eos.ShapeType() )
3602     {
3603     case TopAbs_FACE: {
3604       edge._cosin = 0;
3605       break;
3606     }
3607     case TopAbs_EDGE: {
3608       TopoDS_Edge E    = TopoDS::Edge( eos._shape );
3609       gp_Vec inFaceDir = getFaceDir( F, E, node, helper, normOK );
3610       double angle     = inFaceDir.Angle( edge._normal ); // [0,PI]
3611       edge._cosin      = Cos( angle );
3612       break;
3613     }
3614     case TopAbs_VERTEX: {
3615       if ( fromVonF )
3616       {
3617         getFaceDir( TopoDS::Face( eos._sWOL ), TopoDS::Vertex( eos._shape ),
3618                     node, helper, normOK, &edge._cosin );
3619       }
3620       else if ( eos.SWOLType() != TopAbs_FACE ) // else _cosin is set by getFaceDir()
3621       {
3622         TopoDS_Vertex V  = TopoDS::Vertex( eos._shape );
3623         gp_Vec inFaceDir = getFaceDir( F, V, node, helper, normOK );
3624         double angle     = inFaceDir.Angle( edge._normal ); // [0,PI]
3625         edge._cosin      = Cos( angle );
3626         if ( totalNbFaces > 2 || helper.IsSeamShape( node->getshapeId() ))
3627           for ( int iF = 1; iF < totalNbFaces; ++iF )
3628           {
3629             F = face2Norm[ iF ].first;
3630             inFaceDir = getFaceDir( F, V, node, helper, normOK=true );
3631             if ( normOK ) {
3632               double angle = inFaceDir.Angle( edge._normal );
3633               double cosin = Cos( angle );
3634               if ( Abs( cosin ) > Abs( edge._cosin ))
3635                 edge._cosin = cosin;
3636             }
3637           }
3638       }
3639       break;
3640     }
3641     default:
3642       return error(SMESH_Comment("Invalid shape position of node ")<<node, data._index);
3643     }
3644   }
3645
3646   double normSize = edge._normal.SquareModulus();
3647   if ( normSize < numeric_limits<double>::min() )
3648     return error(SMESH_Comment("Bad normal at node ")<< node->GetID(), data._index );
3649
3650   edge._normal /= sqrt( normSize );
3651
3652   if ( edge.Is( _LayerEdge::MULTI_NORMAL ) && edge._nodes.size() == 2 )
3653   {
3654     getMeshDS()->RemoveFreeNode( edge._nodes.back(), 0, /*fromGroups=*/false );
3655     edge._nodes.resize( 1 );
3656     edge._normal.SetCoord( 0,0,0 );
3657     edge.SetMaxLen( 0 );
3658   }
3659
3660   // Set the rest data
3661   // --------------------
3662
3663   edge.SetCosin( edge._cosin ); // to update edge._lenFactor
3664
3665   if ( onShrinkShape )
3666   {
3667     const SMDS_MeshNode* tgtNode = edge._nodes.back();
3668     if ( SMESHDS_SubMesh* sm = getMeshDS()->MeshElements( data._solid ))
3669       sm->RemoveNode( tgtNode );
3670
3671     // set initial position which is parameters on _sWOL in this case
3672     if ( eos.SWOLType() == TopAbs_EDGE )
3673     {
3674       double u = helper.GetNodeU( TopoDS::Edge( eos._sWOL ), node, 0, &normOK );
3675       edge._pos.push_back( gp_XYZ( u, 0, 0 ));
3676       if ( edge._nodes.size() > 1 )
3677         getMeshDS()->SetNodeOnEdge( tgtNode, TopoDS::Edge( eos._sWOL ), u );
3678     }
3679     else // eos.SWOLType() == TopAbs_FACE
3680     {
3681       gp_XY uv = helper.GetNodeUV( TopoDS::Face( eos._sWOL ), node, 0, &normOK );
3682       edge._pos.push_back( gp_XYZ( uv.X(), uv.Y(), 0));
3683       if ( edge._nodes.size() > 1 )
3684         getMeshDS()->SetNodeOnFace( tgtNode, TopoDS::Face( eos._sWOL ), uv.X(), uv.Y() );
3685     }
3686
3687     if ( edge._nodes.size() > 1 )
3688     {
3689       // check if an angle between a FACE with layers and SWOL is sharp,
3690       // else the edge should not inflate
3691       F.Nullify();
3692       for ( int iF = 0; iF < totalNbFaces  &&  F.IsNull();  ++iF ) // find a FACE with VL
3693         if ( ! helper.IsSubShape( eos._sWOL, face2Norm[iF].first ))
3694           F = face2Norm[iF].first;
3695       if ( !F.IsNull())
3696       {
3697         geomNorm = getFaceNormal( node, F, helper, normOK );
3698         if ( helper.GetSubShapeOri( data._solid, F ) != TopAbs_REVERSED )
3699           geomNorm.Reverse(); // inside the SOLID
3700         if ( geomNorm * edge._normal < -0.001 )
3701         {
3702           getMeshDS()->RemoveFreeNode( tgtNode, 0, /*fromGroups=*/false );
3703           edge._nodes.resize( 1 );
3704         }
3705         else if ( edge._lenFactor > 3 )
3706         {
3707           edge._lenFactor = 2;
3708           edge.Set( _LayerEdge::RISKY_SWOL );
3709         }
3710       }
3711     }
3712   }
3713   else
3714   {
3715     edge._pos.push_back( SMESH_TNodeXYZ( node ));
3716
3717     if ( eos.ShapeType() == TopAbs_FACE )
3718     {
3719       double angle;
3720       for ( size_t i = 0; i < edge._simplices.size(); ++i )
3721       {
3722         edge._simplices[i].IsMinAngleOK( edge._pos.back(), angle );
3723         edge._minAngle = Max( edge._minAngle, angle ); // "angle" is actually cosine
3724       }
3725     }
3726   }
3727
3728   // Set neighbor nodes for a _LayerEdge based on EDGE
3729
3730   if ( eos.ShapeType() == TopAbs_EDGE /*||
3731        ( onShrinkShape && posType == SMDS_TOP_VERTEX && fabs( edge._cosin ) < 1e-10 )*/)
3732   {
3733     edge._2neibors = new _2NearEdges;
3734     // target nodes instead of source ones will be set later
3735   }
3736
3737   return true;
3738 }
3739
3740 //================================================================================
3741 /*!
3742  * \brief Return normal to a FACE at a node
3743  *  \param [in] n - node
3744  *  \param [in] face - FACE
3745  *  \param [in] helper - helper
3746  *  \param [out] isOK - true or false
3747  *  \param [in] shiftInside - to find normal at a position shifted inside the face
3748  *  \return gp_XYZ - normal
3749  */
3750 //================================================================================
3751
3752 gp_XYZ _ViscousBuilder::getFaceNormal(const SMDS_MeshNode* node,
3753                                       const TopoDS_Face&   face,
3754                                       SMESH_MesherHelper&  helper,
3755                                       bool&                isOK,
3756                                       bool                 shiftInside)
3757 {
3758   gp_XY uv;
3759   if ( shiftInside )
3760   {
3761     // get a shifted position
3762     gp_Pnt p = SMESH_TNodeXYZ( node );
3763     gp_XYZ shift( 0,0,0 );
3764     TopoDS_Shape S = helper.GetSubShapeByNode( node, helper.GetMeshDS() );
3765     switch ( S.ShapeType() ) {
3766     case TopAbs_VERTEX:
3767     {
3768       shift = getFaceDir( face, TopoDS::Vertex( S ), node, helper, isOK );
3769       break;
3770     }
3771     case TopAbs_EDGE:
3772     {
3773       shift = getFaceDir( face, TopoDS::Edge( S ), node, helper, isOK );
3774       break;
3775     }
3776     default:
3777       isOK = false;
3778     }
3779     if ( isOK )
3780       shift.Normalize();
3781     p.Translate( shift * 1e-5 );
3782
3783     TopLoc_Location loc;
3784     GeomAPI_ProjectPointOnSurf& projector = helper.GetProjector( face, loc, 1e-7 );
3785
3786     if ( !loc.IsIdentity() ) p.Transform( loc.Transformation().Inverted() );
3787     
3788     projector.Perform( p );
3789     if ( !projector.IsDone() || projector.NbPoints() < 1 )
3790     {
3791       isOK = false;
3792       return p.XYZ();
3793     }
3794     Standard_Real U,V;
3795     projector.LowerDistanceParameters(U,V);
3796     uv.SetCoord( U,V );
3797   }
3798   else
3799   {
3800     uv = helper.GetNodeUV( face, node, 0, &isOK );
3801   }
3802
3803   gp_Dir normal;
3804   isOK = false;
3805
3806   Handle(Geom_Surface) surface = BRep_Tool::Surface( face );
3807
3808   if ( !shiftInside &&
3809        helper.IsDegenShape( node->getshapeId() ) &&
3810        getFaceNormalAtSingularity( uv, face, helper, normal ))
3811   {
3812     isOK = true;
3813     return normal.XYZ();
3814   }
3815
3816   int pointKind = GeomLib::NormEstim( surface, uv, 1e-5, normal );
3817   enum { REGULAR = 0, QUASYSINGULAR, CONICAL, IMPOSSIBLE };
3818
3819   if ( pointKind == IMPOSSIBLE &&
3820        node->GetPosition()->GetDim() == 2 ) // node inside the FACE
3821   {
3822     // probably NormEstim() failed due to a too high tolerance
3823     pointKind = GeomLib::NormEstim( surface, uv, 1e-20, normal );
3824     isOK = ( pointKind < IMPOSSIBLE );
3825   }
3826   if ( pointKind < IMPOSSIBLE )
3827   {
3828     if ( pointKind != REGULAR &&
3829          !shiftInside &&
3830          node->GetPosition()->GetDim() < 2 ) // FACE boundary
3831     {
3832       gp_XYZ normShift = getFaceNormal( node, face, helper, isOK, /*shiftInside=*/true );
3833       if ( normShift * normal.XYZ() < 0. )
3834         normal = normShift;
3835     }
3836     isOK = true;
3837   }
3838
3839   if ( !isOK ) // hard singularity, to call with shiftInside=true ?
3840   {
3841     const TGeomID faceID = helper.GetMeshDS()->ShapeToIndex( face );
3842
3843     SMDS_ElemIteratorPtr fIt = node->GetInverseElementIterator(SMDSAbs_Face);
3844     while ( fIt->more() )
3845     {
3846       const SMDS_MeshElement* f = fIt->next();
3847       if ( f->getshapeId() == faceID )
3848       {
3849         isOK = SMESH_MeshAlgos::FaceNormal( f, (gp_XYZ&) normal.XYZ(), /*normalized=*/true );
3850         if ( isOK )
3851         {
3852           TopoDS_Face ff = face;
3853           ff.Orientation( TopAbs_FORWARD );
3854           if ( helper.IsReversedSubMesh( ff ))
3855             normal.Reverse();
3856           break;
3857         }
3858       }
3859     }
3860   }
3861   return normal.XYZ();
3862 }
3863
3864 //================================================================================
3865 /*!
3866  * \brief Try to get normal at a singularity of a surface basing on it's nature
3867  */
3868 //================================================================================
3869
3870 bool _ViscousBuilder::getFaceNormalAtSingularity( const gp_XY&        uv,
3871                                                   const TopoDS_Face&  face,
3872                                                   SMESH_MesherHelper& helper,
3873                                                   gp_Dir&             normal )
3874 {
3875   BRepAdaptor_Surface surface( face );
3876   gp_Dir axis;
3877   if ( !getRovolutionAxis( surface, axis ))
3878     return false;
3879
3880   double f,l, d, du, dv;
3881   f = surface.FirstUParameter();
3882   l = surface.LastUParameter();
3883   d = ( uv.X() - f ) / ( l - f );
3884   du = ( d < 0.5 ? +1. : -1 ) * 1e-5 * ( l - f );
3885   f = surface.FirstVParameter();
3886   l = surface.LastVParameter();
3887   d = ( uv.Y() - f ) / ( l - f );
3888   dv = ( d < 0.5 ? +1. : -1 ) * 1e-5 * ( l - f );
3889
3890   gp_Dir refDir;
3891   gp_Pnt2d testUV = uv;
3892   enum { REGULAR = 0, QUASYSINGULAR, CONICAL, IMPOSSIBLE };
3893   double tol = 1e-5;
3894   Handle(Geom_Surface) geomsurf = surface.Surface().Surface();
3895   for ( int iLoop = 0; true ; ++iLoop )
3896   {
3897     testUV.SetCoord( testUV.X() + du, testUV.Y() + dv );
3898     if ( GeomLib::NormEstim( geomsurf, testUV, tol, refDir ) == REGULAR )
3899       break;
3900     if ( iLoop > 20 )
3901       return false;
3902     tol /= 10.;
3903   }
3904
3905   if ( axis * refDir < 0. )
3906     axis.Reverse();
3907
3908   normal = axis;
3909
3910   return true;
3911 }
3912
3913 //================================================================================
3914 /*!
3915  * \brief Return a normal at a node weighted with angles taken by faces
3916  */
3917 //================================================================================
3918
3919 gp_XYZ _ViscousBuilder::getWeigthedNormal( const _LayerEdge* edge )
3920 {
3921   const SMDS_MeshNode* n = edge->_nodes[0];
3922
3923   gp_XYZ resNorm(0,0,0);
3924   SMESH_TNodeXYZ p0( n ), pP, pN;
3925   for ( size_t i = 0; i < edge->_simplices.size(); ++i )
3926   {
3927     pP.Set( edge->_simplices[i]._nPrev );
3928     pN.Set( edge->_simplices[i]._nNext );
3929     gp_Vec v0P( p0, pP ), v0N( p0, pN ), vPN( pP, pN ), norm = v0P ^ v0N;
3930     double l0P = v0P.SquareMagnitude();
3931     double l0N = v0N.SquareMagnitude();
3932     double lPN = vPN.SquareMagnitude();
3933     if ( l0P < std::numeric_limits<double>::min() ||
3934          l0N < std::numeric_limits<double>::min() ||
3935          lPN < std::numeric_limits<double>::min() )
3936       continue;
3937     double lNorm = norm.SquareMagnitude();
3938     double  sin2 = lNorm / l0P / l0N;
3939     double angle = ACos(( v0P * v0N ) / Sqrt( l0P ) / Sqrt( l0N ));
3940
3941     double weight = sin2 * angle / lPN;
3942     resNorm += weight * norm.XYZ() / Sqrt( lNorm );
3943   }
3944
3945   return resNorm;
3946 }
3947
3948 //================================================================================
3949 /*!
3950  * \brief Return a normal at a node by getting a common point of offset planes
3951  *        defined by the FACE normals
3952  */
3953 //================================================================================
3954
3955 gp_XYZ _ViscousBuilder::getNormalByOffset( _LayerEdge*                      edge,
3956                                            std::pair< TopoDS_Face, gp_XYZ > f2Normal[],
3957                                            int                              nbFaces,
3958                                            bool                             lastNoOffset)
3959 {
3960   SMESH_TNodeXYZ p0 = edge->_nodes[0];
3961
3962   gp_XYZ resNorm(0,0,0);
3963   TopoDS_Shape V = SMESH_MesherHelper::GetSubShapeByNode( p0._node, getMeshDS() );
3964   if ( V.ShapeType() != TopAbs_VERTEX || nbFaces < 3 )
3965   {
3966     for ( int i = 0; i < nbFaces; ++i )
3967       resNorm += f2Normal[i].second;
3968     return resNorm;
3969   }
3970
3971   // prepare _OffsetPlane's
3972   vector< _OffsetPlane > pln( nbFaces );
3973   for ( int i = 0; i < nbFaces - lastNoOffset; ++i )
3974   {
3975     pln[i]._faceIndex = i;
3976     pln[i]._plane = gp_Pln( p0 + f2Normal[i].second, f2Normal[i].second );
3977   }
3978   if ( lastNoOffset )
3979   {
3980     pln[ nbFaces - 1 ]._faceIndex = nbFaces - 1;
3981     pln[ nbFaces - 1 ]._plane = gp_Pln( p0, f2Normal[ nbFaces - 1 ].second );
3982   }
3983
3984   // intersect neighboring OffsetPlane's
3985   PShapeIteratorPtr edgeIt = SMESH_MesherHelper::GetAncestors( V, *_mesh, TopAbs_EDGE );
3986   while ( const TopoDS_Shape* edge = edgeIt->next() )
3987   {
3988     int f1 = -1, f2 = -1;
3989     for ( int i = 0; i < nbFaces &&  f2 < 0;  ++i )
3990       if ( SMESH_MesherHelper::IsSubShape( *edge, f2Normal[i].first ))
3991         (( f1 < 0 ) ? f1 : f2 ) = i;
3992
3993     if ( f2 >= 0 )
3994       pln[ f1 ].ComputeIntersectionLine( pln[ f2 ], TopoDS::Edge( *edge ), TopoDS::Vertex( V ));
3995   }
3996
3997   // get a common point
3998   gp_XYZ commonPnt( 0, 0, 0 );
3999   int nbPoints = 0;
4000   bool isPointFound;
4001   for ( int i = 0; i < nbFaces; ++i )
4002   {
4003     commonPnt += pln[ i ].GetCommonPoint( isPointFound, TopoDS::Vertex( V ));
4004     nbPoints  += isPointFound;
4005   }
4006   gp_XYZ wgtNorm = getWeigthedNormal( edge );
4007   if ( nbPoints == 0 )
4008     return wgtNorm;
4009
4010   commonPnt /= nbPoints;
4011   resNorm = commonPnt - p0;
4012   if ( lastNoOffset )
4013     return resNorm;
4014
4015   // choose the best among resNorm and wgtNorm
4016   resNorm.Normalize();
4017   wgtNorm.Normalize();
4018   double resMinDot = std::numeric_limits<double>::max();
4019   double wgtMinDot = std::numeric_limits<double>::max();
4020   for ( int i = 0; i < nbFaces - lastNoOffset; ++i )
4021   {
4022     resMinDot = Min( resMinDot, resNorm * f2Normal[i].second );
4023     wgtMinDot = Min( wgtMinDot, wgtNorm * f2Normal[i].second );
4024   }
4025
4026   if ( Max( resMinDot, wgtMinDot ) < theMinSmoothCosin )
4027   {
4028     edge->Set( _LayerEdge::MULTI_NORMAL );
4029   }
4030
4031   return ( resMinDot > wgtMinDot ) ? resNorm : wgtNorm;
4032 }
4033
4034 //================================================================================
4035 /*!
4036  * \brief Compute line of intersection of 2 planes
4037  */
4038 //================================================================================
4039
4040 void _OffsetPlane::ComputeIntersectionLine( _OffsetPlane&        pln,
4041                                             const TopoDS_Edge&   E,
4042                                             const TopoDS_Vertex& V )
4043 {
4044   int iNext = bool( _faceIndexNext[0] >= 0 );
4045   _faceIndexNext[ iNext ] = pln._faceIndex;
4046
4047   gp_XYZ n1 = _plane.Axis().Direction().XYZ();
4048   gp_XYZ n2 = pln._plane.Axis().Direction().XYZ();
4049
4050   gp_XYZ lineDir = n1 ^ n2;
4051
4052   double x = Abs( lineDir.X() );
4053   double y = Abs( lineDir.Y() );
4054   double z = Abs( lineDir.Z() );
4055
4056   int cooMax; // max coordinate
4057   if (x > y) {
4058     if (x > z) cooMax = 1;
4059     else       cooMax = 3;
4060   }
4061   else {
4062     if (y > z) cooMax = 2;
4063     else       cooMax = 3;
4064   }
4065
4066   gp_Pnt linePos;
4067   if ( Abs( lineDir.Coord( cooMax )) < 0.05 )
4068   {
4069     // parallel planes - intersection is an offset of the common EDGE
4070     gp_Pnt p = BRep_Tool::Pnt( V );
4071     linePos  = 0.5 * (( p.XYZ() + n1 ) + ( p.XYZ() + n2 ));
4072     lineDir  = getEdgeDir( E, V );
4073   }
4074   else
4075   {
4076     // the constants in the 2 plane equations
4077     double d1 = - ( _plane.Axis().Direction().XYZ()     * _plane.Location().XYZ() );
4078     double d2 = - ( pln._plane.Axis().Direction().XYZ() * pln._plane.Location().XYZ() );
4079
4080     switch ( cooMax ) {
4081     case 1:
4082       linePos.SetX(  0 );
4083       linePos.SetY(( d2*n1.Z() - d1*n2.Z()) / lineDir.X() );
4084       linePos.SetZ(( d1*n2.Y() - d2*n1.Y()) / lineDir.X() );
4085       break;
4086     case 2:
4087       linePos.SetX(( d1*n2.Z() - d2*n1.Z()) / lineDir.Y() );
4088       linePos.SetY(  0 );
4089       linePos.SetZ(( d2*n1.X() - d1*n2.X()) / lineDir.Y() );
4090       break;
4091     case 3:
4092       linePos.SetX(( d2*n1.Y() - d1*n2.Y()) / lineDir.Z() );
4093       linePos.SetY(( d1*n2.X() - d2*n1.X()) / lineDir.Z() );
4094       linePos.SetZ(  0 );
4095     }
4096   }
4097   gp_Lin& line = _lines[ iNext ];
4098   line.SetDirection( lineDir );
4099   line.SetLocation ( linePos );
4100
4101   _isLineOK[ iNext ] = true;
4102
4103
4104   iNext = bool( pln._faceIndexNext[0] >= 0 );
4105   pln._lines        [ iNext ] = line;
4106   pln._faceIndexNext[ iNext ] = this->_faceIndex;
4107   pln._isLineOK     [ iNext ] = true;
4108 }
4109
4110 //================================================================================
4111 /*!
4112  * \brief Computes intersection point of two _lines
4113  */
4114 //================================================================================
4115
4116 gp_XYZ _OffsetPlane::GetCommonPoint(bool&                 isFound,
4117                                     const TopoDS_Vertex & V) const
4118 {
4119   gp_XYZ p( 0,0,0 );
4120   isFound = false;
4121
4122   if ( NbLines() == 2 )
4123   {
4124     gp_Vec lPerp0 = _lines[0].Direction().XYZ() ^ _plane.Axis().Direction().XYZ();
4125     double  dot01 = lPerp0 * _lines[1].Direction().XYZ();
4126     if ( Abs( dot01 ) > 0.05 )
4127     {
4128       gp_Vec l0l1 = _lines[1].Location().XYZ() - _lines[0].Location().XYZ();
4129       double   u1 = - ( lPerp0 * l0l1 ) / dot01;
4130       p = ( _lines[1].Location().XYZ() + _lines[1].Direction().XYZ() * u1 );
4131       isFound = true;
4132     }
4133     else
4134     {
4135       gp_Pnt  pV ( BRep_Tool::Pnt( V ));
4136       gp_Vec  lv0( _lines[0].Location(), pV    ),  lv1(_lines[1].Location(), pV     );
4137       double dot0( lv0 * _lines[0].Direction() ), dot1( lv1 * _lines[1].Direction() );
4138       p += 0.5 * ( _lines[0].Location().XYZ() + _lines[0].Direction().XYZ() * dot0 );
4139       p += 0.5 * ( _lines[1].Location().XYZ() + _lines[1].Direction().XYZ() * dot1 );
4140       isFound = true;
4141     }
4142   }
4143
4144   return p;
4145 }
4146
4147 //================================================================================
4148 /*!
4149  * \brief Find 2 neighbor nodes of a node on EDGE
4150  */
4151 //================================================================================
4152
4153 bool _ViscousBuilder::findNeiborsOnEdge(const _LayerEdge*     edge,
4154                                         const SMDS_MeshNode*& n1,
4155                                         const SMDS_MeshNode*& n2,
4156                                         _EdgesOnShape&        eos,
4157                                         _SolidData&           data)
4158 {
4159   const SMDS_MeshNode* node = edge->_nodes[0];
4160   const int        shapeInd = eos._shapeID;
4161   SMESHDS_SubMesh*   edgeSM = 0;
4162   if ( eos.ShapeType() == TopAbs_EDGE )
4163   {
4164     edgeSM = eos._subMesh->GetSubMeshDS();
4165     if ( !edgeSM || edgeSM->NbElements() == 0 )
4166       return error(SMESH_Comment("Not meshed EDGE ") << shapeInd, data._index);
4167   }
4168   int iN = 0;
4169   n2 = 0;
4170   SMDS_ElemIteratorPtr eIt = node->GetInverseElementIterator(SMDSAbs_Edge);
4171   while ( eIt->more() && !n2 )
4172   {
4173     const SMDS_MeshElement* e = eIt->next();
4174     const SMDS_MeshNode*   nNeibor = e->GetNode( 0 );
4175     if ( nNeibor == node ) nNeibor = e->GetNode( 1 );
4176     if ( edgeSM )
4177     {
4178       if (!edgeSM->Contains(e)) continue;
4179     }
4180     else
4181     {
4182       TopoDS_Shape s = SMESH_MesherHelper::GetSubShapeByNode( nNeibor, getMeshDS() );
4183       if ( !SMESH_MesherHelper::IsSubShape( s, eos._sWOL )) continue;
4184     }
4185     ( iN++ ? n2 : n1 ) = nNeibor;
4186   }
4187   if ( !n2 )
4188     return error(SMESH_Comment("Wrongly meshed EDGE ") << shapeInd, data._index);
4189   return true;
4190 }
4191
4192 //================================================================================
4193 /*!
4194  * \brief Set _curvature and _2neibors->_plnNorm by 2 neighbor nodes residing the same EDGE
4195  */
4196 //================================================================================
4197
4198 void _LayerEdge::SetDataByNeighbors( const SMDS_MeshNode* n1,
4199                                      const SMDS_MeshNode* n2,
4200                                      const _EdgesOnShape& eos,
4201                                      SMESH_MesherHelper&  helper)
4202 {
4203   if ( eos.ShapeType() != TopAbs_EDGE )
4204     return;
4205   if ( _curvature && Is( SMOOTHED_C1 ))
4206     return;
4207
4208   gp_XYZ  pos = SMESH_TNodeXYZ( _nodes[0] );
4209   gp_XYZ vec1 = pos - SMESH_TNodeXYZ( n1 );
4210   gp_XYZ vec2 = pos - SMESH_TNodeXYZ( n2 );
4211
4212   // Set _curvature
4213
4214   double      sumLen = vec1.Modulus() + vec2.Modulus();
4215   _2neibors->_wgt[0] = 1 - vec1.Modulus() / sumLen;
4216   _2neibors->_wgt[1] = 1 - vec2.Modulus() / sumLen;
4217   double avgNormProj = 0.5 * ( _normal * vec1 + _normal * vec2 );
4218   double      avgLen = 0.5 * ( vec1.Modulus() + vec2.Modulus() );
4219   if ( _curvature ) delete _curvature;
4220   _curvature = _Curvature::New( avgNormProj, avgLen );
4221   // if ( _curvature )
4222   //   debugMsg( _nodes[0]->GetID()
4223   //             << " CURV r,k: " << _curvature->_r<<","<<_curvature->_k
4224   //             << " proj = "<<avgNormProj<< " len = " << avgLen << "| lenDelta(0) = "
4225   //             << _curvature->lenDelta(0) );
4226
4227   // Set _plnNorm
4228
4229   if ( eos._sWOL.IsNull() )
4230   {
4231     TopoDS_Edge  E = TopoDS::Edge( eos._shape );
4232     // if ( SMESH_Algo::isDegenerated( E ))
4233     //   return;
4234     gp_XYZ dirE    = getEdgeDir( E, _nodes[0], helper );
4235     gp_XYZ plnNorm = dirE ^ _normal;
4236     double proj0   = plnNorm * vec1;
4237     double proj1   = plnNorm * vec2;
4238     if ( fabs( proj0 ) > 1e-10 || fabs( proj1 ) > 1e-10 )
4239     {
4240       if ( _2neibors->_plnNorm ) delete _2neibors->_plnNorm;
4241       _2neibors->_plnNorm = new gp_XYZ( plnNorm.Normalized() );
4242     }
4243   }
4244 }
4245
4246 //================================================================================
4247 /*!
4248  * \brief Copy data from a _LayerEdge of other SOLID and based on the same node;
4249  * this and the other _LayerEdge are inflated along a FACE or an EDGE
4250  */
4251 //================================================================================
4252
4253 gp_XYZ _LayerEdge::Copy( _LayerEdge&         other,
4254                          _EdgesOnShape&      eos,
4255                          SMESH_MesherHelper& helper )
4256 {
4257   _nodes     = other._nodes;
4258   _normal    = other._normal;
4259   _len       = 0;
4260   _lenFactor = other._lenFactor;
4261   _cosin     = other._cosin;
4262   _2neibors  = other._2neibors;
4263   _curvature = 0; std::swap( _curvature, other._curvature );
4264   _2neibors  = 0; std::swap( _2neibors,  other._2neibors );
4265
4266   gp_XYZ lastPos( 0,0,0 );
4267   if ( eos.SWOLType() == TopAbs_EDGE )
4268   {
4269     double u = helper.GetNodeU( TopoDS::Edge( eos._sWOL ), _nodes[0] );
4270     _pos.push_back( gp_XYZ( u, 0, 0));
4271
4272     u = helper.GetNodeU( TopoDS::Edge( eos._sWOL ), _nodes.back() );
4273     lastPos.SetX( u );
4274   }
4275   else // TopAbs_FACE
4276   {
4277     gp_XY uv = helper.GetNodeUV( TopoDS::Face( eos._sWOL ), _nodes[0]);
4278     _pos.push_back( gp_XYZ( uv.X(), uv.Y(), 0));
4279
4280     uv = helper.GetNodeUV( TopoDS::Face( eos._sWOL ), _nodes.back() );
4281     lastPos.SetX( uv.X() );
4282     lastPos.SetY( uv.Y() );
4283   }
4284   return lastPos;
4285 }
4286
4287 //================================================================================
4288 /*!
4289  * \brief Set _cosin and _lenFactor
4290  */
4291 //================================================================================
4292
4293 void _LayerEdge::SetCosin( double cosin )
4294 {
4295   _cosin = cosin;
4296   cosin = Abs( _cosin );
4297   //_lenFactor = ( cosin < 1.-1e-12 ) ?  Min( 2., 1./sqrt(1-cosin*cosin )) : 1.0;
4298   _lenFactor = ( cosin < 1.-1e-12 ) ?  1./sqrt(1-cosin*cosin ) : 1.0;
4299 }
4300
4301 //================================================================================
4302 /*!
4303  * \brief Check if another _LayerEdge is a neighbor on EDGE
4304  */
4305 //================================================================================
4306
4307 bool _LayerEdge::IsNeiborOnEdge( const _LayerEdge* edge ) const
4308 {
4309   return (( this->_2neibors && this->_2neibors->include( edge )) ||
4310           ( edge->_2neibors && edge->_2neibors->include( this )));
4311 }
4312
4313 //================================================================================
4314 /*!
4315  * \brief Fills a vector<_Simplex > 
4316  */
4317 //================================================================================
4318
4319 void _Simplex::GetSimplices( const SMDS_MeshNode* node,
4320                              vector<_Simplex>&    simplices,
4321                              const set<TGeomID>&  ingnoreShapes,
4322                              const _SolidData*    dataToCheckOri,
4323                              const bool           toSort)
4324 {
4325   simplices.clear();
4326   SMDS_ElemIteratorPtr fIt = node->GetInverseElementIterator(SMDSAbs_Face);
4327   while ( fIt->more() )
4328   {
4329     const SMDS_MeshElement* f = fIt->next();
4330     const TGeomID    shapeInd = f->getshapeId();
4331     if ( ingnoreShapes.count( shapeInd )) continue;
4332     const int nbNodes = f->NbCornerNodes();
4333     const int  srcInd = f->GetNodeIndex( node );
4334     const SMDS_MeshNode* nPrev = f->GetNode( SMESH_MesherHelper::WrapIndex( srcInd-1, nbNodes ));
4335     const SMDS_MeshNode* nNext = f->GetNode( SMESH_MesherHelper::WrapIndex( srcInd+1, nbNodes ));
4336     const SMDS_MeshNode* nOpp  = f->GetNode( SMESH_MesherHelper::WrapIndex( srcInd+2, nbNodes ));
4337     if ( dataToCheckOri && dataToCheckOri->_reversedFaceIds.count( shapeInd ))
4338       std::swap( nPrev, nNext );
4339     simplices.push_back( _Simplex( nPrev, nNext, ( nbNodes == 3 ? 0 : nOpp )));
4340   }
4341
4342   if ( toSort )
4343     SortSimplices( simplices );
4344 }
4345
4346 //================================================================================
4347 /*!
4348  * \brief Set neighbor simplices side by side
4349  */
4350 //================================================================================
4351
4352 void _Simplex::SortSimplices(vector<_Simplex>& simplices)
4353 {
4354   vector<_Simplex> sortedSimplices( simplices.size() );
4355   sortedSimplices[0] = simplices[0];
4356   size_t nbFound = 0;
4357   for ( size_t i = 1; i < simplices.size(); ++i )
4358   {
4359     for ( size_t j = 1; j < simplices.size(); ++j )
4360       if ( sortedSimplices[i-1]._nNext == simplices[j]._nPrev )
4361       {
4362         sortedSimplices[i] = simplices[j];
4363         nbFound++;
4364         break;
4365       }
4366   }
4367   if ( nbFound == simplices.size() - 1 )
4368     simplices.swap( sortedSimplices );
4369 }
4370
4371 //================================================================================
4372 /*!
4373  * \brief DEBUG. Create groups containing temporary data of _LayerEdge's
4374  */
4375 //================================================================================
4376
4377 void _ViscousBuilder::makeGroupOfLE()
4378 {
4379 #ifdef _DEBUG_
4380   for ( size_t i = 0 ; i < _sdVec.size(); ++i )
4381   {
4382     if ( _sdVec[i]._n2eMap.empty() ) continue;
4383
4384     dumpFunction( SMESH_Comment("make_LayerEdge_") << i );
4385     TNode2Edge::iterator n2e;
4386     for ( n2e = _sdVec[i]._n2eMap.begin(); n2e != _sdVec[i]._n2eMap.end(); ++n2e )
4387     {
4388       _LayerEdge* le = n2e->second;
4389       // for ( size_t iN = 1; iN < le->_nodes.size(); ++iN )
4390       //   dumpCmd(SMESH_Comment("mesh.AddEdge([ ") <<le->_nodes[iN-1]->GetID()
4391       //           << ", " << le->_nodes[iN]->GetID() <<"])");
4392       if ( le ) {
4393         dumpCmd(SMESH_Comment("mesh.AddEdge([ ") <<le->_nodes[0]->GetID()
4394                 << ", " << le->_nodes.back()->GetID() <<"]) # " << le->_flags );
4395       }
4396     }
4397     dumpFunctionEnd();
4398
4399     dumpFunction( SMESH_Comment("makeNormals") << i );
4400     for ( n2e = _sdVec[i]._n2eMap.begin(); n2e != _sdVec[i]._n2eMap.end(); ++n2e )
4401     {
4402       _LayerEdge* edge = n2e->second;
4403       SMESH_TNodeXYZ nXYZ( edge->_nodes[0] );
4404       nXYZ += edge->_normal * _sdVec[i]._stepSize;
4405       dumpCmd(SMESH_Comment("mesh.AddEdge([ ") << edge->_nodes[0]->GetID()
4406               << ", mesh.AddNode( "<< nXYZ.X()<<","<< nXYZ.Y()<<","<< nXYZ.Z()<<")])");
4407     }
4408     dumpFunctionEnd();
4409
4410     dumpFunction( SMESH_Comment("makeTmpFaces_") << i );
4411     dumpCmd( "faceId1 = mesh.NbElements()" );
4412     TopExp_Explorer fExp( _sdVec[i]._solid, TopAbs_FACE );
4413     for ( ; fExp.More(); fExp.Next() )
4414     {
4415       if ( const SMESHDS_SubMesh* sm = _sdVec[i]._proxyMesh->GetProxySubMesh( fExp.Current() ))
4416       {
4417         if ( sm->NbElements() == 0 ) continue;
4418         SMDS_ElemIteratorPtr fIt = sm->GetElements();
4419         while ( fIt->more())
4420         {
4421           const SMDS_MeshElement* e = fIt->next();
4422           SMESH_Comment cmd("mesh.AddFace([");
4423           for ( int j = 0; j < e->NbCornerNodes(); ++j )
4424             cmd << e->GetNode(j)->GetID() << (j+1 < e->NbCornerNodes() ? ",": "])");
4425           dumpCmd( cmd );
4426         }
4427       }
4428     }
4429     dumpCmd( "faceId2 = mesh.NbElements()" );
4430     dumpCmd( SMESH_Comment( "mesh.MakeGroup( 'tmpFaces_" ) << i << "',"
4431              << "SMESH.FACE, SMESH.FT_RangeOfIds,'=',"
4432              << "'%s-%s' % (faceId1+1, faceId2))");
4433     dumpFunctionEnd();
4434   }
4435 #endif
4436 }
4437
4438 //================================================================================
4439 /*!
4440  * \brief Find maximal _LayerEdge length (layer thickness) limited by geometry
4441  */
4442 //================================================================================
4443
4444 void _ViscousBuilder::computeGeomSize( _SolidData& data )
4445 {
4446   data._geomSize = Precision::Infinite();
4447   double intersecDist;
4448   const SMDS_MeshElement* face;
4449   SMESH_MesherHelper helper( *_mesh );
4450
4451   SMESHUtils::Deleter<SMESH_ElementSearcher> searcher
4452     ( SMESH_MeshAlgos::GetElementSearcher( *getMeshDS(),
4453                                            data._proxyMesh->GetFaces( data._solid )));
4454
4455   for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
4456   {
4457     _EdgesOnShape& eos = data._edgesOnShape[ iS ];
4458     if ( eos._edges.empty() )
4459       continue;
4460     // get neighbor faces, intersection with which should not be considered since
4461     // collisions are avoided by means of smoothing
4462     set< TGeomID > neighborFaces;
4463     if ( eos._hyp.ToSmooth() )
4464     {
4465       SMESH_subMeshIteratorPtr subIt =
4466         eos._subMesh->getDependsOnIterator(/*includeSelf=*/eos.ShapeType() != TopAbs_FACE );
4467       while ( subIt->more() )
4468       {
4469         SMESH_subMesh* sm = subIt->next();
4470         PShapeIteratorPtr fIt = helper.GetAncestors( sm->GetSubShape(), *_mesh, TopAbs_FACE );
4471         while ( const TopoDS_Shape* face = fIt->next() )
4472           neighborFaces.insert( getMeshDS()->ShapeToIndex( *face ));
4473       }
4474     }
4475     // find intersections
4476     double thinkness = eos._hyp.GetTotalThickness();
4477     for ( size_t i = 0; i < eos._edges.size(); ++i )
4478     {
4479       if ( eos._edges[i]->Is( _LayerEdge::BLOCKED )) continue;
4480       eos._edges[i]->SetMaxLen( thinkness );
4481       eos._edges[i]->FindIntersection( *searcher, intersecDist, data._epsilon, eos, &face );
4482       if ( intersecDist > 0 && face )
4483       {
4484         data._geomSize = Min( data._geomSize, intersecDist );
4485         if ( !neighborFaces.count( face->getshapeId() ))
4486           eos[i]->SetMaxLen( Min( thinkness, intersecDist / ( face->GetID() < 0 ? 3. : 2. )));
4487       }
4488     }
4489   }
4490
4491   data._maxThickness = 0;
4492   data._minThickness = 1e100;
4493   list< const StdMeshers_ViscousLayers* >::iterator hyp = data._hyps.begin();
4494   for ( ; hyp != data._hyps.end(); ++hyp )
4495   {
4496     data._maxThickness = Max( data._maxThickness, (*hyp)->GetTotalThickness() );
4497     data._minThickness = Min( data._minThickness, (*hyp)->GetTotalThickness() );
4498   }
4499
4500   // Limit inflation step size by geometry size found by intersecting
4501   // normals of _LayerEdge's with mesh faces
4502   if ( data._stepSize > 0.3 * data._geomSize )
4503     limitStepSize( data, 0.3 * data._geomSize );
4504
4505   if ( data._stepSize > data._minThickness )
4506     limitStepSize( data, data._minThickness );
4507
4508
4509   // -------------------------------------------------------------------------
4510   // Detect _LayerEdge which can't intersect with opposite or neighbor layer,
4511   // so no need in detecting intersection at each inflation step
4512   // -------------------------------------------------------------------------
4513
4514   int nbSteps = data._maxThickness / data._stepSize;
4515   if ( nbSteps < 3 || nbSteps * data._n2eMap.size() < 100000 )
4516     return;
4517
4518   vector< const SMDS_MeshElement* > closeFaces;
4519   int nbDetected = 0;
4520
4521   for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
4522   {
4523     _EdgesOnShape& eos = data._edgesOnShape[ iS ];
4524     if ( eos._edges.empty() || eos.ShapeType() != TopAbs_FACE )
4525       continue;
4526
4527     for ( size_t i = 0; i < eos.size(); ++i )
4528     {
4529       SMESH_NodeXYZ p( eos[i]->_nodes[0] );
4530       double radius = data._maxThickness + 2 * eos[i]->_maxLen;
4531       closeFaces.clear();
4532       searcher->GetElementsInSphere( p, radius, SMDSAbs_Face, closeFaces );
4533
4534       bool toIgnore = true;
4535       for ( size_t iF = 0; iF < closeFaces.size()  && toIgnore; ++iF )
4536         if ( !( toIgnore = ( closeFaces[ iF ]->getshapeId() == eos._shapeID ||
4537                              data._ignoreFaceIds.count( closeFaces[ iF ]->getshapeId() ))))
4538         {
4539           // check if a _LayerEdge will inflate in a direction opposite to a direction
4540           // toward a close face
4541           bool allBehind = true;
4542           for ( int iN = 0; iN < closeFaces[ iF ]->NbCornerNodes()  && allBehind; ++iN )
4543           {
4544             SMESH_NodeXYZ pi( closeFaces[ iF ]->GetNode( iN ));
4545             allBehind = (( pi - p ) * eos[i]->_normal < 0.1 * data._stepSize );
4546           }
4547           toIgnore = allBehind;
4548         }
4549
4550
4551       if ( toIgnore ) // no need to detect intersection
4552       {
4553         eos[i]->Set( _LayerEdge::INTERSECTED );
4554         ++nbDetected;
4555       }
4556     }
4557   }
4558
4559   debugMsg( "Nb LE to intersect " << data._n2eMap.size()-nbDetected << ", ignore " << nbDetected );
4560
4561   return;
4562 }
4563
4564 //================================================================================
4565 /*!
4566  * \brief Increase length of _LayerEdge's to reach the required thickness of layers
4567  */
4568 //================================================================================
4569
4570 bool _ViscousBuilder::inflate(_SolidData& data)
4571 {
4572   SMESH_MesherHelper helper( *_mesh );
4573
4574   const double tgtThick = data._maxThickness;
4575
4576   if ( data._stepSize < 1. )
4577     data._epsilon = data._stepSize * 1e-7;
4578
4579   debugMsg( "-- geomSize = " << data._geomSize << ", stepSize = " << data._stepSize );
4580   _pyDump->Pause();
4581
4582   findCollisionEdges( data, helper );
4583
4584   limitMaxLenByCurvature( data, helper );
4585
4586   _pyDump->Resume();
4587
4588   // limit length of _LayerEdge's around MULTI_NORMAL _LayerEdge's
4589   for ( size_t i = 0; i < data._edgesOnShape.size(); ++i )
4590     if ( data._edgesOnShape[i].ShapeType() == TopAbs_VERTEX &&
4591          data._edgesOnShape[i]._edges.size() > 0 &&
4592          data._edgesOnShape[i]._edges[0]->Is( _LayerEdge::MULTI_NORMAL ))
4593     {
4594       data._edgesOnShape[i]._edges[0]->Unset( _LayerEdge::BLOCKED );
4595       data._edgesOnShape[i]._edges[0]->Block( data );
4596     }
4597
4598   const double safeFactor = ( 2*data._maxThickness < data._geomSize ) ? 1 : theThickToIntersection;
4599
4600   double avgThick = 0, curThick = 0, distToIntersection = Precision::Infinite();
4601   int nbSteps = 0, nbRepeats = 0;
4602   while ( avgThick < 0.99 )
4603   {
4604     // new target length
4605     double prevThick = curThick;
4606     curThick += data._stepSize;
4607     if ( curThick > tgtThick )
4608     {
4609       curThick = tgtThick + tgtThick*( 1.-avgThick ) * nbRepeats;
4610       nbRepeats++;
4611     }
4612
4613     double stepSize = curThick - prevThick;
4614     updateNormalsOfSmoothed( data, helper, nbSteps, stepSize ); // to ease smoothing
4615
4616     // Elongate _LayerEdge's
4617     dumpFunction(SMESH_Comment("inflate")<<data._index<<"_step"<<nbSteps); // debug
4618     for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
4619     {
4620       _EdgesOnShape& eos = data._edgesOnShape[iS];
4621       if ( eos._edges.empty() ) continue;
4622
4623       const double shapeCurThick = Min( curThick, eos._hyp.GetTotalThickness() );
4624       for ( size_t i = 0; i < eos._edges.size(); ++i )
4625       {
4626         eos._edges[i]->SetNewLength( shapeCurThick, eos, helper );
4627       }
4628     }
4629     dumpFunctionEnd();
4630
4631     if ( !updateNormals( data, helper, nbSteps, stepSize )) // to avoid collisions
4632       return false;
4633
4634     // Improve and check quality
4635     if ( !smoothAndCheck( data, nbSteps, distToIntersection ))
4636     {
4637       if ( nbSteps > 0 )
4638       {
4639 #ifdef __NOT_INVALIDATE_BAD_SMOOTH
4640         debugMsg("NOT INVALIDATED STEP!");
4641         return error("Smoothing failed", data._index);
4642 #endif
4643         dumpFunction(SMESH_Comment("invalidate")<<data._index<<"_step"<<nbSteps); // debug
4644         for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
4645         {
4646           _EdgesOnShape& eos = data._edgesOnShape[iS];
4647           for ( size_t i = 0; i < eos._edges.size(); ++i )
4648             eos._edges[i]->InvalidateStep( nbSteps+1, eos );
4649         }
4650         dumpFunctionEnd();
4651       }
4652       break; // no more inflating possible
4653     }
4654     nbSteps++;
4655
4656     // Evaluate achieved thickness
4657     avgThick = 0;
4658     int nbActiveEdges = 0;
4659     for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
4660     {
4661       _EdgesOnShape& eos = data._edgesOnShape[iS];
4662       if ( eos._edges.empty() ) continue;
4663
4664       const double shapeTgtThick = eos._hyp.GetTotalThickness();
4665       for ( size_t i = 0; i < eos._edges.size(); ++i )
4666       {
4667         if ( eos._edges[i]->_nodes.size() > 1 )
4668           avgThick    += Min( 1., eos._edges[i]->_len / shapeTgtThick );
4669         else
4670           avgThick    += shapeTgtThick;
4671         nbActiveEdges += ( ! eos._edges[i]->Is( _LayerEdge::BLOCKED ));
4672       }
4673     }
4674     avgThick /= data._n2eMap.size();
4675     debugMsg( "-- Thickness " << curThick << " ("<< avgThick*100 << "%) reached" );
4676
4677 #ifdef BLOCK_INFLATION
4678     if ( nbActiveEdges == 0 )
4679     {
4680       debugMsg( "-- Stop inflation since all _LayerEdge's BLOCKED " );
4681       break;
4682     }
4683 #else
4684     if ( distToIntersection < tgtThick * avgThick * safeFactor && avgThick < 0.9 )
4685     {
4686       debugMsg( "-- Stop inflation since "
4687                 << " distToIntersection( "<<distToIntersection<<" ) < avgThick( "
4688                 << tgtThick * avgThick << " ) * " << safeFactor );
4689       break;
4690     }
4691 #endif
4692
4693     // new step size
4694     limitStepSize( data, 0.25 * distToIntersection );
4695     if ( data._stepSizeNodes[0] )
4696       data._stepSize = data._stepSizeCoeff *
4697         SMESH_TNodeXYZ(data._stepSizeNodes[0]).Distance(data._stepSizeNodes[1]);
4698
4699   } // while ( avgThick < 0.99 )
4700
4701   if ( nbSteps == 0 )
4702     return error("failed at the very first inflation step", data._index);
4703
4704   if ( avgThick < 0.99 )
4705   {
4706     if ( !data._proxyMesh->_warning || data._proxyMesh->_warning->IsOK() )
4707     {
4708       data._proxyMesh->_warning.reset
4709         ( new SMESH_ComputeError (COMPERR_WARNING,
4710                                   SMESH_Comment("Thickness ") << tgtThick <<
4711                                   " of viscous layers not reached,"
4712                                   " average reached thickness is " << avgThick*tgtThick));
4713     }
4714   }
4715
4716   // Restore position of src nodes moved by inflation on _noShrinkShapes
4717   dumpFunction(SMESH_Comment("restoNoShrink_So")<<data._index); // debug
4718   for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
4719   {
4720     _EdgesOnShape& eos = data._edgesOnShape[iS];
4721     if ( !eos._edges.empty() && eos._edges[0]->_nodes.size() == 1 )
4722       for ( size_t i = 0; i < eos._edges.size(); ++i )
4723       {
4724         restoreNoShrink( *eos._edges[ i ] );
4725       }
4726   }
4727   dumpFunctionEnd();
4728
4729   return safeFactor > 0; // == true (avoid warning: unused variable 'safeFactor')
4730 }
4731
4732 //================================================================================
4733 /*!
4734  * \brief Improve quality of layer inner surface and check intersection
4735  */
4736 //================================================================================
4737
4738 bool _ViscousBuilder::smoothAndCheck(_SolidData& data,
4739                                      const int   infStep,
4740                                      double &    distToIntersection)
4741 {
4742   if ( data._nbShapesToSmooth == 0 )
4743     return true; // no shapes needing smoothing
4744
4745   bool moved, improved;
4746   double vol;
4747   vector< _LayerEdge* >    movedEdges, badEdges;
4748   vector< _EdgesOnShape* > eosC1; // C1 continues shapes
4749   vector< bool >           isConcaveFace;
4750
4751   SMESH_MesherHelper helper(*_mesh);
4752   Handle(ShapeAnalysis_Surface) surface;
4753   TopoDS_Face F;
4754
4755   for ( int isFace = 0; isFace < 2; ++isFace ) // smooth on [ EDGEs, FACEs ]
4756   {
4757     const TopAbs_ShapeEnum shapeType = isFace ? TopAbs_FACE : TopAbs_EDGE;
4758
4759     for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
4760     {
4761       _EdgesOnShape& eos = data._edgesOnShape[ iS ];
4762       if ( !eos._toSmooth ||
4763            eos.ShapeType() != shapeType ||
4764            eos._edges.empty() )
4765         continue;
4766
4767       // already smoothed?
4768       // bool toSmooth = ( eos._edges[ 0 ]->NbSteps() >= infStep+1 );
4769       // if ( !toSmooth ) continue;
4770
4771       if ( !eos._hyp.ToSmooth() )
4772       {
4773         // smooth disabled by the user; check validy only
4774         if ( !isFace ) continue;
4775         badEdges.clear();
4776         for ( size_t i = 0; i < eos._edges.size(); ++i )
4777         {
4778           _LayerEdge* edge = eos._edges[i];
4779           for ( size_t iF = 0; iF < edge->_simplices.size(); ++iF )
4780             if ( !edge->_simplices[iF].IsForward( edge->_nodes[0], edge->_pos.back(), vol ))
4781             {
4782               // debugMsg( "-- Stop inflation. Bad simplex ("
4783               //           << " "<< edge->_nodes[0]->GetID()
4784               //           << " "<< edge->_nodes.back()->GetID()
4785               //           << " "<< edge->_simplices[iF]._nPrev->GetID()
4786               //           << " "<< edge->_simplices[iF]._nNext->GetID() << " ) ");
4787               // return false;
4788               badEdges.push_back( edge );
4789             }
4790         }
4791         if ( !badEdges.empty() )
4792         {
4793           eosC1.resize(1);
4794           eosC1[0] = &eos;
4795           int nbBad = invalidateBadSmooth( data, helper, badEdges, eosC1, infStep );
4796           if ( nbBad > 0 )
4797             return false;
4798         }
4799         continue; // goto the next EDGE or FACE
4800       }
4801
4802       // prepare data
4803       if ( eos.SWOLType() == TopAbs_FACE )
4804       {
4805         if ( !F.IsSame( eos._sWOL )) {
4806           F = TopoDS::Face( eos._sWOL );
4807           helper.SetSubShape( F );
4808           surface = helper.GetSurface( F );
4809         }
4810       }
4811       else
4812       {
4813         F.Nullify(); surface.Nullify();
4814       }
4815       const TGeomID sInd = eos._shapeID;
4816
4817       // perform smoothing
4818
4819       if ( eos.ShapeType() == TopAbs_EDGE )
4820       {
4821         dumpFunction(SMESH_Comment("smooth")<<data._index << "_Ed"<<sInd <<"_InfStep"<<infStep);
4822
4823         if ( !eos._edgeSmoother->Perform( data, surface, F, helper ))
4824         {
4825           // smooth on EDGE's (normally we should not get here)
4826           int step = 0;
4827           do {
4828             moved = false;
4829             for ( size_t i = 0; i < eos._edges.size(); ++i )
4830             {
4831               moved |= eos._edges[i]->SmoothOnEdge( surface, F, helper );
4832             }
4833             dumpCmd( SMESH_Comment("# end step ")<<step);
4834           }
4835           while ( moved && step++ < 5 );
4836         }
4837         dumpFunctionEnd();
4838       }
4839
4840       else // smooth on FACE
4841       {
4842         eosC1.clear();
4843         eosC1.push_back( & eos );
4844         eosC1.insert( eosC1.end(), eos._eosC1.begin(), eos._eosC1.end() );
4845
4846         movedEdges.clear();
4847         isConcaveFace.resize( eosC1.size() );
4848         for ( size_t iEOS = 0; iEOS < eosC1.size(); ++iEOS )
4849         {
4850           isConcaveFace[ iEOS ] = data._concaveFaces.count( eosC1[ iEOS ]->_shapeID  );
4851           vector< _LayerEdge* > & edges = eosC1[ iEOS ]->_edges;
4852           for ( size_t i = 0; i < edges.size(); ++i )
4853             if ( edges[i]->Is( _LayerEdge::MOVED ) ||
4854                  edges[i]->Is( _LayerEdge::NEAR_BOUNDARY ))
4855               movedEdges.push_back( edges[i] );
4856
4857           makeOffsetSurface( *eosC1[ iEOS ], helper );
4858         }
4859
4860         int step = 0, stepLimit = 5, nbBad = 0;
4861         while (( ++step <= stepLimit ) || improved )
4862         {
4863           dumpFunction(SMESH_Comment("smooth")<<data._index<<"_Fa"<<sInd
4864                        <<"_InfStep"<<infStep<<"_"<<step); // debug
4865           int oldBadNb = nbBad;
4866           badEdges.clear();
4867
4868 #ifdef INCREMENTAL_SMOOTH
4869           bool findBest = false; // ( step == stepLimit );
4870           for ( size_t i = 0; i < movedEdges.size(); ++i )
4871           {
4872             movedEdges[i]->Unset( _LayerEdge::SMOOTHED );
4873             if ( movedEdges[i]->Smooth( step, findBest, movedEdges ) > 0 )
4874               badEdges.push_back( movedEdges[i] );
4875           }
4876 #else
4877           bool findBest = ( step == stepLimit || isConcaveFace[ iEOS ]);
4878           for ( size_t iEOS = 0; iEOS < eosC1.size(); ++iEOS )
4879           {
4880             vector< _LayerEdge* > & edges = eosC1[ iEOS ]->_edges;
4881             for ( size_t i = 0; i < edges.size(); ++i )
4882             {
4883               edges[i]->Unset( _LayerEdge::SMOOTHED );
4884               if ( edges[i]->Smooth( step, findBest, false ) > 0 )
4885                 badEdges.push_back( eos._edges[i] );
4886             }
4887           }
4888 #endif
4889           nbBad = badEdges.size();
4890
4891           if ( nbBad > 0 )
4892             debugMsg(SMESH_Comment("nbBad = ") << nbBad );
4893
4894           if ( !badEdges.empty() && step >= stepLimit / 2 )
4895           {
4896             if ( badEdges[0]->Is( _LayerEdge::ON_CONCAVE_FACE ))
4897               stepLimit = 9;
4898
4899             // resolve hard smoothing situation around concave VERTEXes
4900             for ( size_t iEOS = 0; iEOS < eosC1.size(); ++iEOS )
4901             {
4902               vector< _EdgesOnShape* > & eosCoVe = eosC1[ iEOS ]->_eosConcaVer;
4903               for ( size_t i = 0; i < eosCoVe.size(); ++i )
4904                 eosCoVe[i]->_edges[0]->MoveNearConcaVer( eosCoVe[i], eosC1[ iEOS ],
4905                                                          step, badEdges );
4906             }
4907             // look for the best smooth of _LayerEdge's neighboring badEdges
4908             nbBad = 0;
4909             for ( size_t i = 0; i < badEdges.size(); ++i )
4910             {
4911               _LayerEdge* ledge = badEdges[i];
4912               for ( size_t iN = 0; iN < ledge->_neibors.size(); ++iN )
4913               {
4914                 ledge->_neibors[iN]->Unset( _LayerEdge::SMOOTHED );
4915                 nbBad += ledge->_neibors[iN]->Smooth( step, true, /*findBest=*/true );
4916               }
4917               ledge->Unset( _LayerEdge::SMOOTHED );
4918               nbBad += ledge->Smooth( step, true, /*findBest=*/true );
4919             }
4920             debugMsg(SMESH_Comment("nbBad = ") << nbBad );
4921           }
4922
4923           if ( nbBad == oldBadNb  &&
4924                nbBad > 0 &&
4925                step < stepLimit ) // smooth w/o check of validity
4926           {
4927             dumpFunctionEnd();
4928             dumpFunction(SMESH_Comment("smoothWoCheck")<<data._index<<"_Fa"<<sInd
4929                          <<"_InfStep"<<infStep<<"_"<<step); // debug
4930             for ( size_t i = 0; i < movedEdges.size(); ++i )
4931             {
4932               movedEdges[i]->SmoothWoCheck();
4933             }
4934             if ( stepLimit < 9 )
4935               stepLimit++;
4936           }
4937
4938           improved = ( nbBad < oldBadNb );
4939
4940           dumpFunctionEnd();
4941
4942           if (( step % 3 == 1 ) || ( nbBad > 0 && step >= stepLimit / 2 ))
4943             for ( size_t iEOS = 0; iEOS < eosC1.size(); ++iEOS )
4944             {
4945               putOnOffsetSurface( *eosC1[ iEOS ], infStep, eosC1, step, /*moveAll=*/step == 1 );
4946             }
4947
4948         } // smoothing steps
4949
4950         // project -- to prevent intersections or fix bad simplices
4951         for ( size_t iEOS = 0; iEOS < eosC1.size(); ++iEOS )
4952         {
4953           if ( ! eosC1[ iEOS ]->_eosConcaVer.empty() || nbBad > 0 )
4954             putOnOffsetSurface( *eosC1[ iEOS ], infStep, eosC1 );
4955         }
4956
4957         //if ( !badEdges.empty() )
4958         {
4959           badEdges.clear();
4960           for ( size_t iEOS = 0; iEOS < eosC1.size(); ++iEOS )
4961           {
4962             for ( size_t i = 0; i < eosC1[ iEOS ]->_edges.size(); ++i )
4963             {
4964               if ( !eosC1[ iEOS ]->_sWOL.IsNull() ) continue;
4965
4966               _LayerEdge* edge = eosC1[ iEOS ]->_edges[i];
4967               edge->CheckNeiborsOnBoundary( & badEdges );
4968               if (( nbBad > 0 ) ||
4969                   ( edge->Is( _LayerEdge::BLOCKED ) && edge->Is( _LayerEdge::NEAR_BOUNDARY )))
4970               {
4971                 SMESH_TNodeXYZ tgtXYZ = edge->_nodes.back();
4972                 gp_XYZ        prevXYZ = edge->PrevCheckPos();
4973                 for ( size_t j = 0; j < edge->_simplices.size(); ++j )
4974                   if ( !edge->_simplices[j].IsForward( &prevXYZ, &tgtXYZ, vol ))
4975                   {
4976                     debugMsg("Bad simplex ( " << edge->_nodes[0]->GetID()
4977                              << " "<< tgtXYZ._node->GetID()
4978                              << " "<< edge->_simplices[j]._nPrev->GetID()
4979                              << " "<< edge->_simplices[j]._nNext->GetID() << " )" );
4980                     badEdges.push_back( edge );
4981                     break;
4982                   }
4983               }
4984             }
4985           }
4986
4987           // try to fix bad simplices by removing the last inflation step of some _LayerEdge's
4988           nbBad = invalidateBadSmooth( data, helper, badEdges, eosC1, infStep );
4989
4990           if ( nbBad > 0 )
4991             return false;
4992         }
4993
4994       } // // smooth on FACE's
4995     } // loop on shapes
4996   } // smooth on [ EDGEs, FACEs ]
4997
4998   // Check orientation of simplices of _LayerEdge's on EDGEs and VERTEXes
4999   eosC1.resize(1);
5000   for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
5001   {
5002     _EdgesOnShape& eos = data._edgesOnShape[ iS ];
5003     if ( eos.ShapeType() == TopAbs_FACE ||
5004          eos._edges.empty() ||
5005          !eos._sWOL.IsNull() )
5006       continue;
5007
5008     badEdges.clear();
5009     for ( size_t i = 0; i < eos._edges.size(); ++i )
5010     {
5011       _LayerEdge*      edge = eos._edges[i];
5012       if ( edge->_nodes.size() < 2 ) continue;
5013       SMESH_TNodeXYZ tgtXYZ = edge->_nodes.back();
5014       //SMESH_TNodeXYZ prevXYZ = edge->_nodes[0];
5015       gp_XYZ        prevXYZ = edge->PrevCheckPos( &eos );
5016       //const gp_XYZ& prevXYZ = edge->PrevPos();
5017       for ( size_t j = 0; j < edge->_simplices.size(); ++j )
5018         if ( !edge->_simplices[j].IsForward( &prevXYZ, &tgtXYZ, vol ))
5019         {
5020           debugMsg("Bad simplex on bnd ( " << edge->_nodes[0]->GetID()
5021                    << " "<< tgtXYZ._node->GetID()
5022                    << " "<< edge->_simplices[j]._nPrev->GetID()
5023                    << " "<< edge->_simplices[j]._nNext->GetID() << " )" );
5024           badEdges.push_back( edge );
5025           break;
5026         }
5027     }
5028
5029     // try to fix bad simplices by removing the last inflation step of some _LayerEdge's
5030     eosC1[0] = &eos;
5031     int nbBad = invalidateBadSmooth( data, helper, badEdges, eosC1, infStep );
5032     if ( nbBad > 0 )
5033       return false;
5034   }
5035
5036
5037   // Check if the last segments of _LayerEdge intersects 2D elements;
5038   // checked elements are either temporary faces or faces on surfaces w/o the layers
5039
5040   SMESHUtils::Deleter<SMESH_ElementSearcher> searcher
5041     ( SMESH_MeshAlgos::GetElementSearcher( *getMeshDS(),
5042                                            data._proxyMesh->GetFaces( data._solid )) );
5043
5044 #ifdef BLOCK_INFLATION
5045   const bool toBlockInfaltion = true;
5046 #else
5047   const bool toBlockInfaltion = false;
5048 #endif
5049   distToIntersection = Precision::Infinite();
5050   double dist;
5051   const SMDS_MeshElement* intFace = 0;
5052   const SMDS_MeshElement* closestFace = 0;
5053   _LayerEdge* le = 0;
5054   bool is1stBlocked = true; // dbg
5055   for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
5056   {
5057     _EdgesOnShape& eos = data._edgesOnShape[ iS ];
5058     if ( eos._edges.empty() || !eos._sWOL.IsNull() )
5059       continue;
5060     for ( size_t i = 0; i < eos._edges.size(); ++i )
5061     {
5062       if ( eos._edges[i]->Is( _LayerEdge::INTERSECTED ) ||
5063            eos._edges[i]->Is( _LayerEdge::MULTI_NORMAL ))
5064         continue;
5065       if ( eos._edges[i]->FindIntersection( *searcher, dist, data._epsilon, eos, &intFace ))
5066       {
5067         return false;
5068         // commented due to "Illegal hash-positionPosition" error in NETGEN
5069         // on Debian60 on viscous_layers_01/B2 case
5070         // Collision; try to deflate _LayerEdge's causing it
5071         // badEdges.clear();
5072         // badEdges.push_back( eos._edges[i] );
5073         // eosC1[0] = & eos;
5074         // int nbBad = invalidateBadSmooth( data, helper, badEdges, eosC1, infStep );
5075         // if ( nbBad > 0 )
5076         //   return false;
5077
5078         // badEdges.clear();
5079         // if ( _EdgesOnShape* eof = data.GetShapeEdges( intFace->getshapeId() ))
5080         // {
5081         //   if ( const _TmpMeshFace* f = dynamic_cast< const _TmpMeshFace*>( intFace ))
5082         //   {
5083         //     const SMDS_MeshElement* srcFace =
5084         //       eof->_subMesh->GetSubMeshDS()->GetElement( f->getIdInShape() );
5085         //     SMDS_ElemIteratorPtr nIt = srcFace->nodesIterator();
5086         //     while ( nIt->more() )
5087         //     {
5088         //       const SMDS_MeshNode* srcNode = static_cast<const SMDS_MeshNode*>( nIt->next() );
5089         //       TNode2Edge::iterator n2e = data._n2eMap.find( srcNode );
5090         //       if ( n2e != data._n2eMap.end() )
5091         //         badEdges.push_back( n2e->second );
5092         //     }
5093         //     eosC1[0] = eof;
5094         //     nbBad = invalidateBadSmooth( data, helper, badEdges, eosC1, infStep );
5095         //     if ( nbBad > 0 )
5096         //       return false;
5097         //   }
5098         // }
5099         // if ( eos._edges[i]->FindIntersection( *searcher, dist, data._epsilon, eos, &intFace ))
5100         //   return false;
5101         // else
5102         //   continue;
5103       }
5104       if ( !intFace )
5105       {
5106         SMESH_Comment msg("Invalid? normal at node "); msg << eos._edges[i]->_nodes[0]->GetID();
5107         debugMsg( msg );
5108         continue;
5109       }
5110
5111       const bool isShorterDist = ( distToIntersection > dist );
5112       if ( toBlockInfaltion || isShorterDist )
5113       {
5114         // ignore intersection of a _LayerEdge based on a _ConvexFace with a face
5115         // lying on this _ConvexFace
5116         if ( _ConvexFace* convFace = data.GetConvexFace( intFace->getshapeId() ))
5117           if ( convFace->_isTooCurved && convFace->_subIdToEOS.count ( eos._shapeID ))
5118             continue;
5119
5120         // ignore intersection of a _LayerEdge based on a FACE with an element on this FACE
5121         // ( avoid limiting the thickness on the case of issue 22576)
5122         if ( intFace->getshapeId() == eos._shapeID  )
5123           continue;
5124
5125         // ignore intersection with intFace of an adjacent FACE
5126         if ( dist > 0.1 * eos._edges[i]->_len )
5127         {
5128           bool toIgnore = false;
5129           if (  eos._toSmooth )
5130           {
5131             const TopoDS_Shape& S = getMeshDS()->IndexToShape( intFace->getshapeId() );
5132             if ( !S.IsNull() && S.ShapeType() == TopAbs_FACE )
5133             {
5134               TopExp_Explorer sub( eos._shape,
5135                                    eos.ShapeType() == TopAbs_FACE ? TopAbs_EDGE : TopAbs_VERTEX );
5136               for ( ; !toIgnore && sub.More(); sub.Next() )
5137                 // is adjacent - has a common EDGE or VERTEX
5138                 toIgnore = ( helper.IsSubShape( sub.Current(), S ));
5139
5140               if ( toIgnore ) // check angle between normals
5141               {
5142                 gp_XYZ normal;
5143                 if ( SMESH_MeshAlgos::FaceNormal( intFace, normal, /*normalized=*/true ))
5144                   toIgnore  = ( normal * eos._edges[i]->_normal > -0.5 );
5145               }
5146             }
5147           }
5148           if ( !toIgnore ) // check if the edge is a neighbor of intFace
5149           {
5150             for ( size_t iN = 0; !toIgnore &&  iN < eos._edges[i]->_neibors.size(); ++iN )
5151             {
5152               int nInd = intFace->GetNodeIndex( eos._edges[i]->_neibors[ iN ]->_nodes.back() );
5153               toIgnore = ( nInd >= 0 );
5154             }
5155           }
5156           if ( toIgnore )
5157             continue;
5158         }
5159
5160         // intersection not ignored
5161
5162         if ( toBlockInfaltion &&
5163              dist < ( eos._edges[i]->_len * theThickToIntersection ))
5164         {
5165           if ( is1stBlocked ) { is1stBlocked = false; // debug
5166             dumpFunction(SMESH_Comment("blockIntersected") <<data._index<<"_InfStep"<<infStep);
5167           }
5168           eos._edges[i]->Set( _LayerEdge::INTERSECTED ); // not to intersect
5169           eos._edges[i]->Block( data );                  // not to inflate
5170
5171           //if ( _EdgesOnShape* eof = data.GetShapeEdges( intFace->getshapeId() ))
5172           {
5173             // block _LayerEdge's, on top of which intFace is
5174             if ( const _TmpMeshFace* f = dynamic_cast< const _TmpMeshFace*>( intFace ))
5175             {
5176               const SMDS_MeshElement* srcFace = f->_srcFace;
5177               SMDS_ElemIteratorPtr        nIt = srcFace->nodesIterator();
5178               while ( nIt->more() )
5179               {
5180                 const SMDS_MeshNode* srcNode = static_cast<const SMDS_MeshNode*>( nIt->next() );
5181                 TNode2Edge::iterator n2e = data._n2eMap.find( srcNode );
5182                 if ( n2e != data._n2eMap.end() )
5183                   n2e->second->Block( data );
5184               }
5185             }
5186           }
5187         }
5188
5189         if ( isShorterDist )
5190         {
5191           distToIntersection = dist;
5192           le = eos._edges[i];
5193           closestFace = intFace;
5194         }
5195
5196       } // if ( toBlockInfaltion || isShorterDist )
5197     } // loop on eos._edges
5198   } // loop on data._edgesOnShape
5199
5200   if ( !is1stBlocked )
5201     dumpFunctionEnd();
5202
5203   if ( closestFace && le )
5204   {
5205 #ifdef __myDEBUG
5206     SMDS_MeshElement::iterator nIt = closestFace->begin_nodes();
5207     cout << "#Shortest distance: _LayerEdge nodes: tgt " << le->_nodes.back()->GetID()
5208          << " src " << le->_nodes[0]->GetID()<< ", intersection with face ("
5209          << (*nIt++)->GetID()<<" "<< (*nIt++)->GetID()<<" "<< (*nIt++)->GetID()
5210          << ") distance = " << distToIntersection<< endl;
5211 #endif
5212   }
5213
5214   return true;
5215 }
5216
5217 //================================================================================
5218 /*!
5219  * \brief try to fix bad simplices by removing the last inflation step of some _LayerEdge's
5220  *  \param [in,out] badSmooEdges - _LayerEdge's to fix
5221  *  \return int - resulting nb of bad _LayerEdge's
5222  */
5223 //================================================================================
5224
5225 int _ViscousBuilder::invalidateBadSmooth( _SolidData&               data,
5226                                           SMESH_MesherHelper&       helper,
5227                                           vector< _LayerEdge* >&    badSmooEdges,
5228                                           vector< _EdgesOnShape* >& eosC1,
5229                                           const int                 infStep )
5230 {
5231   if ( badSmooEdges.empty() || infStep == 0 ) return 0;
5232
5233   dumpFunction(SMESH_Comment("invalidateBadSmooth")<<"_S"<<eosC1[0]->_shapeID<<"_InfStep"<<infStep);
5234
5235   enum {
5236     INVALIDATED   = _LayerEdge::UNUSED_FLAG,
5237     TO_INVALIDATE = _LayerEdge::UNUSED_FLAG * 2,
5238     ADDED         = _LayerEdge::UNUSED_FLAG * 4
5239   };
5240   data.UnmarkEdges( TO_INVALIDATE & INVALIDATED & ADDED );
5241
5242   double vol;
5243   bool haveInvalidated = true;
5244   while ( haveInvalidated )
5245   {
5246     haveInvalidated = false;
5247     for ( size_t i = 0; i < badSmooEdges.size(); ++i )
5248     {
5249       _LayerEdge*   edge = badSmooEdges[i];
5250       _EdgesOnShape* eos = data.GetShapeEdges( edge );
5251       edge->Set( ADDED );
5252       bool invalidated = false;
5253       if ( edge->Is( TO_INVALIDATE ) && edge->NbSteps() > 1 )
5254       {
5255         edge->InvalidateStep( edge->NbSteps(), *eos, /*restoreLength=*/true );
5256         edge->Block( data );
5257         edge->Set( INVALIDATED );
5258         edge->Unset( TO_INVALIDATE );
5259         invalidated = true;
5260         haveInvalidated = true;
5261       }
5262
5263       // look for _LayerEdge's of bad _simplices
5264       int nbBad = 0;
5265       SMESH_TNodeXYZ tgtXYZ  = edge->_nodes.back();
5266       gp_XYZ        prevXYZ1 = edge->PrevCheckPos( eos );
5267       //const gp_XYZ& prevXYZ2 = edge->PrevPos();
5268       for ( size_t j = 0; j < edge->_simplices.size(); ++j )
5269       {
5270         if (( edge->_simplices[j].IsForward( &prevXYZ1, &tgtXYZ, vol ))/* &&
5271             ( &prevXYZ1 == &prevXYZ2 || edge->_simplices[j].IsForward( &prevXYZ2, &tgtXYZ, vol ))*/)
5272           continue;
5273
5274         bool isBad = true;
5275         _LayerEdge* ee[2] = { 0,0 };
5276         for ( size_t iN = 0; iN < edge->_neibors.size() &&   !ee[1]  ; ++iN )
5277           if ( edge->_simplices[j].Includes( edge->_neibors[iN]->_nodes.back() ))
5278             ee[ ee[0] != 0 ] = edge->_neibors[iN];
5279
5280         int maxNbSteps = Max( ee[0]->NbSteps(), ee[1]->NbSteps() );
5281         while ( maxNbSteps > edge->NbSteps() && isBad )
5282         {
5283           --maxNbSteps;
5284           for ( int iE = 0; iE < 2; ++iE )
5285           {
5286             if ( ee[ iE ]->NbSteps() > maxNbSteps &&
5287                  ee[ iE ]->NbSteps() > 1 )
5288             {
5289               _EdgesOnShape* eos = data.GetShapeEdges( ee[ iE ] );
5290               ee[ iE ]->InvalidateStep( ee[ iE ]->NbSteps(), *eos, /*restoreLength=*/true );
5291               ee[ iE ]->Block( data );
5292               ee[ iE ]->Set( INVALIDATED );
5293               haveInvalidated = true;
5294             }
5295           }
5296           if (( edge->_simplices[j].IsForward( &prevXYZ1, &tgtXYZ, vol )) /*&&
5297               ( &prevXYZ1 == &prevXYZ2 || edge->_simplices[j].IsForward( &prevXYZ2, &tgtXYZ, vol ))*/)
5298             isBad = false;
5299         }
5300         nbBad += isBad;
5301         if ( !ee[0]->Is( ADDED )) badSmooEdges.push_back( ee[0] );
5302         if ( !ee[1]->Is( ADDED )) badSmooEdges.push_back( ee[1] );
5303         ee[0]->Set( ADDED );
5304         ee[1]->Set( ADDED );
5305         if ( isBad )
5306         {
5307           ee[0]->Set( TO_INVALIDATE );
5308           ee[1]->Set( TO_INVALIDATE );
5309         }
5310       }
5311
5312       if ( !invalidated &&  nbBad > 0  &&  edge->NbSteps() > 1 )
5313       {
5314         _EdgesOnShape* eos = data.GetShapeEdges( edge );
5315         edge->InvalidateStep( edge->NbSteps(), *eos, /*restoreLength=*/true );
5316         edge->Block( data );
5317         edge->Set( INVALIDATED );
5318         edge->Unset( TO_INVALIDATE );
5319         haveInvalidated = true;
5320       }
5321     } // loop on badSmooEdges
5322   } // while ( haveInvalidated )
5323
5324   // re-smooth on analytical EDGEs
5325   for ( size_t i = 0; i < badSmooEdges.size(); ++i )
5326   {
5327     _LayerEdge* edge = badSmooEdges[i];
5328     if ( !edge->Is( INVALIDATED )) continue;
5329
5330     _EdgesOnShape* eos = data.GetShapeEdges( edge );
5331     if ( eos->ShapeType() == TopAbs_VERTEX )
5332     {
5333       PShapeIteratorPtr eIt = helper.GetAncestors( eos->_shape, *_mesh, TopAbs_EDGE );
5334       while ( const TopoDS_Shape* e = eIt->next() )
5335         if ( _EdgesOnShape* eoe = data.GetShapeEdges( *e ))
5336           if ( eoe->_edgeSmoother && eoe->_edgeSmoother->isAnalytic() )
5337           {
5338             // TopoDS_Face F; Handle(ShapeAnalysis_Surface) surface;
5339             // if ( eoe->SWOLType() == TopAbs_FACE ) {
5340             //   F       = TopoDS::Face( eoe->_sWOL );
5341             //   surface = helper.GetSurface( F );
5342             // }
5343             // eoe->_edgeSmoother->Perform( data, surface, F, helper );
5344             eoe->_edgeSmoother->_anaCurve.Nullify();
5345           }
5346     }
5347   }
5348
5349
5350   // check result of invalidation
5351
5352   int nbBad = 0;
5353   for ( size_t iEOS = 0; iEOS < eosC1.size(); ++iEOS )
5354   {
5355     for ( size_t i = 0; i < eosC1[ iEOS ]->_edges.size(); ++i )
5356     {
5357       if ( !eosC1[ iEOS ]->_sWOL.IsNull() ) continue;
5358       _LayerEdge*      edge = eosC1[ iEOS ]->_edges[i];
5359       SMESH_TNodeXYZ tgtXYZ = edge->_nodes.back();
5360       gp_XYZ        prevXYZ = edge->PrevCheckPos( eosC1[ iEOS ]);
5361       for ( size_t j = 0; j < edge->_simplices.size(); ++j )
5362         if ( !edge->_simplices[j].IsForward( &prevXYZ, &tgtXYZ, vol ))
5363         {
5364           ++nbBad;
5365           debugMsg("Bad simplex remains ( " << edge->_nodes[0]->GetID()
5366                    << " "<< tgtXYZ._node->GetID()
5367                    << " "<< edge->_simplices[j]._nPrev->GetID()
5368                    << " "<< edge->_simplices[j]._nNext->GetID() << " )" );
5369         }
5370     }
5371   }
5372   dumpFunctionEnd();
5373
5374   return nbBad;
5375 }
5376
5377 //================================================================================
5378 /*!
5379  * \brief Create an offset surface
5380  */
5381 //================================================================================
5382
5383 void _ViscousBuilder::makeOffsetSurface( _EdgesOnShape& eos, SMESH_MesherHelper& helper )
5384 {
5385   if ( eos._offsetSurf.IsNull() ||
5386        eos._edgeForOffset == 0 ||
5387        eos._edgeForOffset->Is( _LayerEdge::BLOCKED ))
5388     return;
5389
5390   Handle(ShapeAnalysis_Surface) baseSurface = helper.GetSurface( TopoDS::Face( eos._shape ));
5391
5392   // find offset
5393   gp_Pnt   tgtP = SMESH_TNodeXYZ( eos._edgeForOffset->_nodes.back() );
5394   /*gp_Pnt2d uv=*/baseSurface->ValueOfUV( tgtP, Precision::Confusion() );
5395   double offset = baseSurface->Gap();
5396
5397   eos._offsetSurf.Nullify();
5398
5399   try
5400   {
5401     BRepOffsetAPI_MakeOffsetShape offsetMaker;
5402     offsetMaker.PerformByJoin( eos._shape, -offset, Precision::Confusion() );
5403     if ( !offsetMaker.IsDone() ) return;
5404
5405     TopExp_Explorer fExp( offsetMaker.Shape(), TopAbs_FACE );
5406     if ( !fExp.More() ) return;
5407
5408     TopoDS_Face F = TopoDS::Face( fExp.Current() );
5409     Handle(Geom_Surface) surf = BRep_Tool::Surface( F );
5410     if ( surf.IsNull() ) return;
5411
5412     eos._offsetSurf = new ShapeAnalysis_Surface( surf );
5413   }
5414   catch ( Standard_Failure )
5415   {
5416   }
5417 }
5418
5419 //================================================================================
5420 /*!
5421  * \brief Put nodes of a curved FACE to its offset surface
5422  */
5423 //================================================================================
5424
5425 void _ViscousBuilder::putOnOffsetSurface( _EdgesOnShape&            eos,
5426                                           int                       infStep,
5427                                           vector< _EdgesOnShape* >& eosC1,
5428                                           int                       smooStep,
5429                                           int                       moveAll )
5430 {
5431   _EdgesOnShape * eof = & eos;
5432   if ( eos.ShapeType() != TopAbs_FACE ) // eos is a boundary of C1 FACE, look for the FACE eos
5433   {
5434     eof = 0;
5435     for ( size_t i = 0; i < eosC1.size() && !eof; ++i )
5436     {
5437       if ( eosC1[i]->_offsetSurf.IsNull() ||
5438            eosC1[i]->ShapeType() != TopAbs_FACE ||
5439            eosC1[i]->_edgeForOffset == 0 ||
5440            eosC1[i]->_edgeForOffset->Is( _LayerEdge::BLOCKED ))
5441         continue;
5442       if ( SMESH_MesherHelper::IsSubShape( eos._shape, eosC1[i]->_shape ))
5443         eof = eosC1[i];
5444     }
5445   }
5446   if ( !eof ||
5447        eof->_offsetSurf.IsNull() ||
5448        eof->ShapeType() != TopAbs_FACE ||
5449        eof->_edgeForOffset == 0 ||
5450        eof->_edgeForOffset->Is( _LayerEdge::BLOCKED ))
5451     return;
5452
5453   double preci = BRep_Tool::Tolerance( TopoDS::Face( eof->_shape )), vol;
5454   for ( size_t i = 0; i < eos._edges.size(); ++i )
5455   {
5456     _LayerEdge* edge = eos._edges[i];
5457     edge->Unset( _LayerEdge::MARKED );
5458     if ( edge->Is( _LayerEdge::BLOCKED ) || !edge->_curvature )
5459       continue;
5460     if ( moveAll == _LayerEdge::UPD_NORMAL_CONV )
5461     {
5462       if ( !edge->Is( _LayerEdge::UPD_NORMAL_CONV ))
5463         continue;
5464     }
5465     else if ( !moveAll && !edge->Is( _LayerEdge::MOVED ))
5466       continue;
5467
5468     int nbBlockedAround = 0;
5469     for ( size_t iN = 0; iN < edge->_neibors.size(); ++iN )
5470       nbBlockedAround += edge->_neibors[iN]->Is( _LayerEdge::BLOCKED );
5471     if ( nbBlockedAround > 1 )
5472       continue;
5473
5474     gp_Pnt tgtP = SMESH_TNodeXYZ( edge->_nodes.back() );
5475     gp_Pnt2d uv = eof->_offsetSurf->NextValueOfUV( edge->_curvature->_uv, tgtP, preci );
5476     if ( eof->_offsetSurf->Gap() > edge->_len ) continue; // NextValueOfUV() bug
5477     edge->_curvature->_uv = uv;
5478     if ( eof->_offsetSurf->Gap() < 10 * preci ) continue; // same pos
5479
5480     gp_XYZ  newP = eof->_offsetSurf->Value( uv ).XYZ();
5481     gp_XYZ prevP = edge->PrevCheckPos();
5482     bool      ok = true;
5483     if ( !moveAll )
5484       for ( size_t iS = 0; iS < edge->_simplices.size() && ok; ++iS )
5485       {
5486         ok = edge->_simplices[iS].IsForward( &prevP, &newP, vol );
5487       }
5488     if ( ok )
5489     {
5490       SMDS_MeshNode* n = const_cast< SMDS_MeshNode* >( edge->_nodes.back() );
5491       n->setXYZ( newP.X(), newP.Y(), newP.Z());
5492       edge->_pos.back() = newP;
5493
5494       edge->Set( _LayerEdge::MARKED );
5495       if ( moveAll == _LayerEdge::UPD_NORMAL_CONV )
5496       {
5497         edge->_normal = ( newP - prevP ).Normalized();
5498       }
5499     }
5500   }
5501
5502
5503
5504 #ifdef _DEBUG_
5505   // dumpMove() for debug
5506   size_t i = 0;
5507   for ( ; i < eos._edges.size(); ++i )
5508     if ( eos._edges[i]->Is( _LayerEdge::MARKED ))
5509       break;
5510   if ( i < eos._edges.size() )
5511   {
5512     dumpFunction(SMESH_Comment("putOnOffsetSurface_S") << eos._shapeID
5513                  << "_InfStep" << infStep << "_" << smooStep );
5514     for ( ; i < eos._edges.size(); ++i )
5515     {
5516       if ( eos._edges[i]->Is( _LayerEdge::MARKED ))
5517         dumpMove( eos._edges[i]->_nodes.back() );
5518     }
5519     dumpFunctionEnd();
5520   }
5521 #endif
5522
5523   _ConvexFace* cnvFace;
5524   if ( moveAll != _LayerEdge::UPD_NORMAL_CONV &&
5525        eos.ShapeType() == TopAbs_FACE &&
5526        (cnvFace = eos.GetData().GetConvexFace( eos._shapeID )) &&
5527        !cnvFace->_normalsFixedOnBorders )
5528   {
5529     // put on the surface nodes built on FACE boundaries
5530     SMESH_subMeshIteratorPtr smIt = eos._subMesh->getDependsOnIterator(/*includeSelf=*/false);
5531     while ( smIt->more() )
5532     {
5533       SMESH_subMesh* sm = smIt->next();
5534       _EdgesOnShape* subEOS = eos.GetData().GetShapeEdges( sm->GetId() );
5535       if ( !subEOS->_sWOL.IsNull() ) continue;
5536       if ( std::find( eosC1.begin(), eosC1.end(), subEOS ) != eosC1.end() ) continue;
5537
5538       putOnOffsetSurface( *subEOS, infStep, eosC1, smooStep, _LayerEdge::UPD_NORMAL_CONV );
5539     }
5540     cnvFace->_normalsFixedOnBorders = true;
5541   }
5542 }
5543
5544 //================================================================================
5545 /*!
5546  * \brief Return a curve of the EDGE to be used for smoothing and arrange
5547  *        _LayerEdge's to be in a consequent order
5548  */
5549 //================================================================================
5550
5551 Handle(Geom_Curve) _Smoother1D::CurveForSmooth( const TopoDS_Edge&  E,
5552                                                 _EdgesOnShape&      eos,
5553                                                 SMESH_MesherHelper& helper)
5554 {
5555   SMESHDS_SubMesh* smDS = eos._subMesh->GetSubMeshDS();
5556
5557   TopLoc_Location loc; double f,l;
5558
5559   Handle(Geom_Line)   line;
5560   Handle(Geom_Circle) circle;
5561   bool isLine, isCirc;
5562   if ( eos._sWOL.IsNull() ) /////////////////////////////////////////// 3D case
5563   {
5564     // check if the EDGE is a line
5565     Handle(Geom_Curve) curve = BRep_Tool::Curve( E, f, l);
5566     if ( curve->IsKind( STANDARD_TYPE( Geom_TrimmedCurve )))
5567       curve = Handle(Geom_TrimmedCurve)::DownCast( curve )->BasisCurve();
5568
5569     line   = Handle(Geom_Line)::DownCast( curve );
5570     circle = Handle(Geom_Circle)::DownCast( curve );
5571     isLine = (!line.IsNull());
5572     isCirc = (!circle.IsNull());
5573
5574     if ( !isLine && !isCirc ) // Check if the EDGE is close to a line
5575     {
5576       isLine = SMESH_Algo::IsStraight( E );
5577
5578       if ( isLine )
5579         line = new Geom_Line( gp::OX() ); // only type does matter
5580     }
5581     if ( !isLine && !isCirc && eos._edges.size() > 2) // Check if the EDGE is close to a circle
5582     {
5583       // TODO
5584     }
5585   }
5586   else //////////////////////////////////////////////////////////////////////// 2D case
5587   {
5588     if ( !eos._isRegularSWOL ) // 23190
5589       return NULL;
5590
5591     const TopoDS_Face& F = TopoDS::Face( eos._sWOL );
5592
5593     // check if the EDGE is a line
5594     Handle(Geom2d_Curve) curve = BRep_Tool::CurveOnSurface( E, F, f, l );
5595     if ( curve->IsKind( STANDARD_TYPE( Geom2d_TrimmedCurve )))
5596       curve = Handle(Geom2d_TrimmedCurve)::DownCast( curve )->BasisCurve();
5597
5598     Handle(Geom2d_Line)   line2d   = Handle(Geom2d_Line)::DownCast( curve );
5599     Handle(Geom2d_Circle) circle2d = Handle(Geom2d_Circle)::DownCast( curve );
5600     isLine = (!line2d.IsNull());
5601     isCirc = (!circle2d.IsNull());
5602
5603     if ( !isLine && !isCirc ) // Check if the EDGE is close to a line
5604     {
5605       Bnd_B2d bndBox;
5606       SMDS_NodeIteratorPtr nIt = smDS->GetNodes();
5607       while ( nIt->more() )
5608         bndBox.Add( helper.GetNodeUV( F, nIt->next() ));
5609       gp_XY size = bndBox.CornerMax() - bndBox.CornerMin();
5610
5611       const double lineTol = 1e-2 * sqrt( bndBox.SquareExtent() );
5612       for ( int i = 0; i < 2 && !isLine; ++i )
5613         isLine = ( size.Coord( i+1 ) <= lineTol );
5614     }
5615     if ( !isLine && !isCirc && eos._edges.size() > 2 ) // Check if the EDGE is close to a circle
5616     {
5617       // TODO
5618     }
5619     if ( isLine )
5620     {
5621       line = new Geom_Line( gp::OX() ); // only type does matter
5622     }
5623     else if ( isCirc )
5624     {
5625       gp_Pnt2d p = circle2d->Location();
5626       gp_Ax2 ax( gp_Pnt( p.X(), p.Y(), 0), gp::DX());
5627       circle = new Geom_Circle( ax, 1.); // only center position does matter
5628     }
5629   }
5630
5631   if ( isLine )
5632     return line;
5633   if ( isCirc )
5634     return circle;
5635
5636   return Handle(Geom_Curve)();
5637 }
5638
5639 //================================================================================
5640 /*!
5641  * \brief Smooth edges on EDGE
5642  */
5643 //================================================================================
5644
5645 bool _Smoother1D::Perform(_SolidData&                    data,
5646                           Handle(ShapeAnalysis_Surface)& surface,
5647                           const TopoDS_Face&             F,
5648                           SMESH_MesherHelper&            helper )
5649 {
5650   if ( _leParams.empty() || ( !isAnalytic() && _offPoints.empty() ))
5651     prepare( data );
5652
5653   findEdgesToSmooth();
5654   if ( isAnalytic() )
5655     return smoothAnalyticEdge( data, surface, F, helper );
5656   else
5657     return smoothComplexEdge ( data, surface, F, helper );
5658 }
5659
5660 //================================================================================
5661 /*!
5662  * \brief Find edges to smooth
5663  */
5664 //================================================================================
5665
5666 void _Smoother1D::findEdgesToSmooth()
5667 {
5668   _LayerEdge* leOnV[2] = { getLEdgeOnV(0), getLEdgeOnV(1) };
5669   for ( int iEnd = 0; iEnd < 2; ++iEnd )
5670     if ( leOnV[iEnd]->Is( _LayerEdge::NORMAL_UPDATED ))
5671       _leOnV[iEnd]._cosin = Abs( _edgeDir[iEnd].Normalized() * leOnV[iEnd]->_normal );
5672
5673   _eToSmooth[0].first = _eToSmooth[0].second = 0;
5674
5675   for ( size_t i = 0; i < _eos.size(); ++i )
5676   {
5677     if ( !_eos[i]->Is( _LayerEdge::TO_SMOOTH ))
5678     {
5679       if ( needSmoothing( _leOnV[0]._cosin,
5680                           _eos[i]->_len * leOnV[0]->_lenFactor, _curveLen * _leParams[i] ) ||
5681            isToSmooth( i )
5682            )
5683         _eos[i]->Set( _LayerEdge::TO_SMOOTH );
5684       else
5685         break;
5686     }
5687     _eToSmooth[0].second = i+1;
5688   }
5689
5690   _eToSmooth[1].first = _eToSmooth[1].second = _eos.size();
5691
5692   for ( int i = _eos.size() - 1; i >= _eToSmooth[0].second; --i )
5693   {
5694     if ( !_eos[i]->Is( _LayerEdge::TO_SMOOTH ))
5695     {
5696       if ( needSmoothing( _leOnV[1]._cosin,
5697                           _eos[i]->_len * leOnV[1]->_lenFactor, _curveLen * ( 1.-_leParams[i] )) ||
5698            isToSmooth( i ))
5699         _eos[i]->Set( _LayerEdge::TO_SMOOTH );
5700       else
5701         break;
5702     }
5703     _eToSmooth[1].first = i;
5704   }
5705 }
5706
5707 //================================================================================
5708 /*!
5709  * \brief Check if iE-th _LayerEdge needs smoothing
5710  */
5711 //================================================================================
5712
5713 bool _Smoother1D::isToSmooth( int iE )
5714 {
5715   SMESH_NodeXYZ pi( _eos[iE]->_nodes[0] );
5716   SMESH_NodeXYZ p0( _eos[iE]->_2neibors->srcNode(0) );
5717   SMESH_NodeXYZ p1( _eos[iE]->_2neibors->srcNode(1) );
5718   gp_XYZ       seg0 = pi - p0;
5719   gp_XYZ       seg1 = p1 - pi;
5720   gp_XYZ    tangent =  seg0 + seg1;
5721   double tangentLen = tangent.Modulus();
5722   double  segMinLen = Min( seg0.Modulus(), seg1.Modulus() );
5723   if ( tangentLen < std::numeric_limits<double>::min() )
5724     return false;
5725   tangent /= tangentLen;
5726
5727   for ( size_t i = 0; i < _eos[iE]->_neibors.size(); ++i )
5728   {
5729     _LayerEdge* ne = _eos[iE]->_neibors[i];
5730     if ( !ne->Is( _LayerEdge::TO_SMOOTH ) ||
5731          ne->_nodes.size() < 2 ||
5732          ne->_nodes[0]->GetPosition()->GetDim() != 2 )
5733       continue;
5734     gp_XYZ edgeVec = SMESH_NodeXYZ( ne->_nodes.back() ) - SMESH_NodeXYZ( ne->_nodes[0] );
5735     double    proj = edgeVec * tangent;
5736     if ( needSmoothing( 1., proj, segMinLen ))
5737       return true;
5738   }
5739   return false;
5740 }
5741
5742 //================================================================================
5743 /*!
5744  * \brief smooth _LayerEdge's on a staight EDGE or circular EDGE
5745  */
5746 //================================================================================
5747
5748 bool _Smoother1D::smoothAnalyticEdge( _SolidData&                    data,
5749                                       Handle(ShapeAnalysis_Surface)& surface,
5750                                       const TopoDS_Face&             F,
5751                                       SMESH_MesherHelper&            helper)
5752 {
5753   if ( !isAnalytic() ) return false;
5754
5755   size_t iFrom = 0, iTo = _eos._edges.size();
5756
5757   if ( _anaCurve->IsKind( STANDARD_TYPE( Geom_Line )))
5758   {
5759     if ( F.IsNull() ) // 3D
5760     {
5761       SMESH_TNodeXYZ pSrc0( _eos._edges[iFrom]->_2neibors->srcNode(0) );
5762       SMESH_TNodeXYZ pSrc1( _eos._edges[iTo-1]->_2neibors->srcNode(1) );
5763       //const   gp_XYZ lineDir = pSrc1 - pSrc0;
5764       //_LayerEdge* vLE0 = getLEdgeOnV( 0 );
5765       //_LayerEdge* vLE1 = getLEdgeOnV( 1 );
5766       // bool shiftOnly = ( vLE0->Is( _LayerEdge::NORMAL_UPDATED ) ||
5767       //                    vLE0->Is( _LayerEdge::BLOCKED ) ||
5768       //                    vLE1->Is( _LayerEdge::NORMAL_UPDATED ) ||
5769       //                    vLE1->Is( _LayerEdge::BLOCKED ));
5770       for ( int iEnd = 0; iEnd < 2; ++iEnd )
5771       {
5772         iFrom = _eToSmooth[ iEnd ].first, iTo = _eToSmooth[ iEnd ].second;
5773         if ( iFrom >= iTo ) continue;
5774         SMESH_TNodeXYZ p0( _eos[iFrom]->_2neibors->tgtNode(0) );
5775         SMESH_TNodeXYZ p1( _eos[iTo-1]->_2neibors->tgtNode(1) );
5776         double param0 = ( iFrom == 0 ) ? 0. : _leParams[ iFrom-1 ];
5777         double param1 = _leParams[ iTo ];
5778         for ( size_t i = iFrom; i < iTo; ++i )
5779         {
5780           _LayerEdge*       edge = _eos[i];
5781           SMDS_MeshNode* tgtNode = const_cast<SMDS_MeshNode*>( edge->_nodes.back() );
5782           double           param = ( _leParams[i] - param0 ) / ( param1 - param0 );
5783           gp_XYZ          newPos = p0 * ( 1. - param ) + p1 * param;
5784
5785           // if ( shiftOnly || edge->Is( _LayerEdge::NORMAL_UPDATED ))
5786           // {
5787           //   gp_XYZ curPos = SMESH_TNodeXYZ ( tgtNode );
5788           //   double  shift = ( lineDir * ( newPos - pSrc0 ) -
5789           //                     lineDir * ( curPos - pSrc0 ));
5790           //   newPos = curPos + lineDir * shift / lineDir.SquareModulus();
5791           // }
5792           if ( edge->Is( _LayerEdge::BLOCKED ))
5793           {
5794             SMESH_TNodeXYZ pSrc( edge->_nodes[0] );
5795             double curThick = pSrc.SquareDistance( tgtNode );
5796             double newThink = ( pSrc - newPos ).SquareModulus();
5797             if ( newThink > curThick )
5798               continue;
5799           }
5800           edge->_pos.back() = newPos;
5801           tgtNode->setXYZ( newPos.X(), newPos.Y(), newPos.Z() );
5802           dumpMove( tgtNode );
5803         }
5804       }
5805     }
5806     else // 2D
5807     {
5808       _LayerEdge* eV0 = getLEdgeOnV( 0 );
5809       _LayerEdge* eV1 = getLEdgeOnV( 1 );
5810       gp_XY      uvV0 = eV0->LastUV( F, *data.GetShapeEdges( eV0 ));
5811       gp_XY      uvV1 = eV1->LastUV( F, *data.GetShapeEdges( eV1 ));
5812       if ( eV0->_nodes.back() == eV1->_nodes.back() ) // closed edge
5813       {
5814         int iPeriodic = helper.GetPeriodicIndex();
5815         if ( iPeriodic == 1 || iPeriodic == 2 )
5816         {
5817           uvV1.SetCoord( iPeriodic, helper.GetOtherParam( uvV1.Coord( iPeriodic )));
5818           if ( uvV0.Coord( iPeriodic ) > uvV1.Coord( iPeriodic ))
5819             std::swap( uvV0, uvV1 );
5820         }
5821       }
5822       for ( int iEnd = 0; iEnd < 2; ++iEnd )
5823       {
5824         iFrom = _eToSmooth[ iEnd ].first, iTo = _eToSmooth[ iEnd ].second;
5825         if ( iFrom >= iTo ) continue;
5826         _LayerEdge* e0 = _eos[iFrom]->_2neibors->_edges[0];
5827         _LayerEdge* e1 = _eos[iTo-1]->_2neibors->_edges[1];
5828         gp_XY      uv0 = ( e0 == eV0 ) ? uvV0 : e0->LastUV( F, _eos );
5829         gp_XY      uv1 = ( e1 == eV1 ) ? uvV1 : e1->LastUV( F, _eos );
5830         double  param0 = ( iFrom == 0 ) ? 0. : _leParams[ iFrom-1 ];
5831         double  param1 = _leParams[ iTo ];
5832         gp_XY  rangeUV = uv1 - uv0;
5833         for ( size_t i = iFrom; i < iTo; ++i )
5834         {
5835           if ( _eos[i]->Is( _LayerEdge::BLOCKED )) continue;
5836           double param = ( _leParams[i] - param0 ) / ( param1 - param0 );
5837           gp_XY newUV = uv0 + param * rangeUV;
5838
5839           gp_Pnt newPos = surface->Value( newUV.X(), newUV.Y() );
5840           SMDS_MeshNode* tgtNode = const_cast<SMDS_MeshNode*>( _eos[i]->_nodes.back() );
5841           tgtNode->setXYZ( newPos.X(), newPos.Y(), newPos.Z() );
5842           dumpMove( tgtNode );
5843
5844           SMDS_FacePositionPtr pos = tgtNode->GetPosition();
5845           pos->SetUParameter( newUV.X() );
5846           pos->SetVParameter( newUV.Y() );
5847
5848           gp_XYZ newUV0( newUV.X(), newUV.Y(), 0 );
5849
5850           if ( !_eos[i]->Is( _LayerEdge::SMOOTHED ))
5851           {
5852             _eos[i]->Set( _LayerEdge::SMOOTHED ); // to check in refine() (IPAL54237)
5853             if ( _eos[i]->_pos.size() > 2 )
5854             {
5855               // modify previous positions to make _LayerEdge less sharply bent
5856               vector<gp_XYZ>& uvVec = _eos[i]->_pos;
5857               const gp_XYZ  uvShift = newUV0 - uvVec.back();
5858               const double     len2 = ( uvVec.back() - uvVec[ 0 ] ).SquareModulus();
5859               int iPrev = uvVec.size() - 2;
5860               while ( iPrev > 0 )
5861               {
5862                 double r = ( uvVec[ iPrev ] - uvVec[0] ).SquareModulus() / len2;
5863                 uvVec[ iPrev ] += uvShift * r;
5864                 --iPrev;
5865               }
5866             }
5867           }
5868           _eos[i]->_pos.back() = newUV0;
5869         }
5870       }
5871     }
5872     return true;
5873   }
5874
5875   if ( _anaCurve->IsKind( STANDARD_TYPE( Geom_Circle )))
5876   {
5877     Handle(Geom_Circle) circle = Handle(Geom_Circle)::DownCast( _anaCurve );
5878     gp_Pnt center3D = circle->Location();
5879
5880     if ( F.IsNull() ) // 3D
5881     {
5882       if ( getLEdgeOnV( 0 )->_nodes.back() == getLEdgeOnV( 1 )->_nodes.back() )
5883         return true; // closed EDGE - nothing to do
5884
5885       // circle is a real curve of EDGE
5886       gp_Circ circ = circle->Circ();
5887
5888       // new center is shifted along its axis
5889       const gp_Dir& axis = circ.Axis().Direction();
5890       _LayerEdge*     e0 = getLEdgeOnV(0);
5891       _LayerEdge*     e1 = getLEdgeOnV(1);
5892       SMESH_TNodeXYZ  p0 = e0->_nodes.back();
5893       SMESH_TNodeXYZ  p1 = e1->_nodes.back();
5894       double      shift1 = axis.XYZ() * ( p0 - center3D.XYZ() );
5895       double      shift2 = axis.XYZ() * ( p1 - center3D.XYZ() );
5896       gp_Pnt   newCenter = center3D.XYZ() + axis.XYZ() * 0.5 * ( shift1 + shift2 );
5897
5898       double newRadius = 0.5 * ( newCenter.Distance( p0 ) + newCenter.Distance( p1 ));
5899
5900       gp_Ax2  newAxis( newCenter, axis, gp_Vec( newCenter, p0 ));
5901       gp_Circ newCirc( newAxis, newRadius );
5902       gp_Vec  vecC1  ( newCenter, p1 );
5903
5904       double uLast = newAxis.XDirection().AngleWithRef( vecC1, newAxis.Direction() ); // -PI - +PI
5905       if ( uLast < 0 )
5906         uLast += 2 * M_PI;
5907       
5908       for ( size_t i = 0; i < _eos.size(); ++i )
5909       {
5910         if ( _eos[i]->Is( _LayerEdge::BLOCKED )) continue;
5911         //if ( !_eos[i]->Is( _LayerEdge::TO_SMOOTH )) continue;
5912         double u = uLast * _leParams[i];
5913         gp_Pnt p = ElCLib::Value( u, newCirc );
5914         _eos._edges[i]->_pos.back() = p.XYZ();
5915
5916         SMDS_MeshNode* tgtNode = const_cast<SMDS_MeshNode*>( _eos._edges[i]->_nodes.back() );
5917         tgtNode->setXYZ( p.X(), p.Y(), p.Z() );
5918         dumpMove( tgtNode );
5919       }
5920       return true;
5921     }
5922     else // 2D
5923     {
5924       const gp_XY center( center3D.X(), center3D.Y() );
5925
5926       _LayerEdge* e0 = getLEdgeOnV(0);
5927       _LayerEdge* eM = _eos._edges[ 0 ];
5928       _LayerEdge* e1 = getLEdgeOnV(1);
5929       gp_XY      uv0 = e0->LastUV( F, *data.GetShapeEdges( e0 ) );
5930       gp_XY      uvM = eM->LastUV( F, *data.GetShapeEdges( eM ) );
5931       gp_XY      uv1 = e1->LastUV( F, *data.GetShapeEdges( e1 ) );
5932       gp_Vec2d vec0( center, uv0 );
5933       gp_Vec2d vecM( center, uvM );
5934       gp_Vec2d vec1( center, uv1 );
5935       double uLast = vec0.Angle( vec1 ); // -PI - +PI
5936       double uMidl = vec0.Angle( vecM );
5937       if ( uLast * uMidl <= 0. )
5938         uLast += ( uMidl > 0 ? +2. : -2. ) * M_PI;
5939       const double radius = 0.5 * ( vec0.Magnitude() + vec1.Magnitude() );
5940
5941       gp_Ax2d   axis( center, vec0 );
5942       gp_Circ2d circ( axis, radius );
5943       for ( size_t i = 0; i < _eos.size(); ++i )
5944       {
5945         if ( _eos[i]->Is( _LayerEdge::BLOCKED )) continue;
5946         //if ( !_eos[i]->Is( _LayerEdge::TO_SMOOTH )) continue;
5947         double    newU = uLast * _leParams[i];
5948         gp_Pnt2d newUV = ElCLib::Value( newU, circ );
5949         _eos._edges[i]->_pos.back().SetCoord( newUV.X(), newUV.Y(), 0 );
5950
5951         gp_Pnt newPos = surface->Value( newUV.X(), newUV.Y() );
5952         SMDS_MeshNode* tgtNode = const_cast<SMDS_MeshNode*>( _eos._edges[i]->_nodes.back() );
5953         tgtNode->setXYZ( newPos.X(), newPos.Y(), newPos.Z() );
5954         dumpMove( tgtNode );
5955
5956         SMDS_FacePositionPtr pos = tgtNode->GetPosition();
5957         pos->SetUParameter( newUV.X() );
5958         pos->SetVParameter( newUV.Y() );
5959
5960         _eos[i]->Set( _LayerEdge::SMOOTHED ); // to check in refine() (IPAL54237)
5961       }
5962     }
5963     return true;
5964   }
5965
5966   return false;
5967 }
5968
5969 //================================================================================
5970 /*!
5971  * \brief smooth _LayerEdge's on a an EDGE
5972  */
5973 //================================================================================
5974
5975 bool _Smoother1D::smoothComplexEdge( _SolidData&                    data,
5976                                      Handle(ShapeAnalysis_Surface)& surface,
5977                                      const TopoDS_Face&             F,
5978                                      SMESH_MesherHelper&            helper)
5979 {
5980   if ( _offPoints.empty() )
5981     return false;
5982
5983   // ----------------------------------------------
5984   // move _offPoints along normals of _LayerEdge's
5985   // ----------------------------------------------
5986
5987   _LayerEdge* e[2] = { getLEdgeOnV(0), getLEdgeOnV(1) };
5988   if ( e[0]->Is( _LayerEdge::NORMAL_UPDATED ))
5989     _leOnV[0]._normal = getNormalNormal( e[0]->_normal, _edgeDir[0] );
5990   if ( e[1]->Is( _LayerEdge::NORMAL_UPDATED )) 
5991     _leOnV[1]._normal = getNormalNormal( e[1]->_normal, _edgeDir[1] );
5992   _leOnV[0]._len = e[0]->_len;
5993   _leOnV[1]._len = e[1]->_len;
5994   for ( size_t i = 0; i < _offPoints.size(); i++ )
5995   {
5996     _LayerEdge*  e0 = _offPoints[i]._2edges._edges[0];
5997     _LayerEdge*  e1 = _offPoints[i]._2edges._edges[1];
5998     const double w0 = _offPoints[i]._2edges._wgt[0];
5999     const double w1 = _offPoints[i]._2edges._wgt[1];
6000     gp_XYZ  avgNorm = ( e0->_normal    * w0 + e1->_normal    * w1 ).Normalized();
6001     double  avgLen  = ( e0->_len       * w0 + e1->_len       * w1 );
6002     double  avgFact = ( e0->_lenFactor * w0 + e1->_lenFactor * w1 );
6003     if ( e0->Is( _LayerEdge::NORMAL_UPDATED ) ||
6004          e1->Is( _LayerEdge::NORMAL_UPDATED ))
6005       avgNorm = getNormalNormal( avgNorm, _offPoints[i]._edgeDir );
6006
6007     _offPoints[i]._xyz += avgNorm * ( avgLen - _offPoints[i]._len ) * avgFact;
6008     _offPoints[i]._len  = avgLen;
6009   }
6010
6011   double fTol = 0;
6012   if ( !surface.IsNull() ) // project _offPoints to the FACE
6013   {
6014     fTol = 100 * BRep_Tool::Tolerance( F );
6015     //const double segLen = _offPoints[0].Distance( _offPoints[1] );
6016
6017     gp_Pnt2d uv = surface->ValueOfUV( _offPoints[0]._xyz, fTol );
6018     //if ( surface->Gap() < 0.5 * segLen )
6019       _offPoints[0]._xyz = surface->Value( uv ).XYZ();
6020
6021     for ( size_t i = 1; i < _offPoints.size(); ++i )
6022     {
6023       uv = surface->NextValueOfUV( uv, _offPoints[i]._xyz, fTol );
6024       //if ( surface->Gap() < 0.5 * segLen )
6025         _offPoints[i]._xyz = surface->Value( uv ).XYZ();
6026     }
6027   }
6028
6029   // -----------------------------------------------------------------
6030   // project tgt nodes of extreme _LayerEdge's to the offset segments
6031   // -----------------------------------------------------------------
6032
6033   const int updatedOrBlocked = _LayerEdge::NORMAL_UPDATED | _LayerEdge::BLOCKED;
6034   if ( e[0]->Is( updatedOrBlocked )) _iSeg[0] = 0;
6035   if ( e[1]->Is( updatedOrBlocked )) _iSeg[1] = _offPoints.size()-2;
6036
6037   gp_Pnt pExtreme[2], pProj[2];
6038   bool isProjected[2];
6039   for ( int is2nd = 0; is2nd < 2; ++is2nd )
6040   {
6041     pExtreme[ is2nd ] = SMESH_TNodeXYZ( e[is2nd]->_nodes.back() );
6042     int  i = _iSeg[ is2nd ];
6043     int di = is2nd ? -1 : +1;
6044     bool & projected = isProjected[ is2nd ];
6045     projected = false;
6046     double uOnSeg, distMin = Precision::Infinite(), dist, distPrev = 0;
6047     int nbWorse = 0;
6048     do {
6049       gp_Vec v0p( _offPoints[i]._xyz, pExtreme[ is2nd ]    );
6050       gp_Vec v01( _offPoints[i]._xyz, _offPoints[i+1]._xyz );
6051       uOnSeg     = ( v0p * v01 ) / v01.SquareMagnitude();  // param [0,1] along v01
6052       projected  = ( Abs( uOnSeg - 0.5 ) <= 0.5 );
6053       dist       =  pExtreme[ is2nd ].SquareDistance( _offPoints[ i + ( uOnSeg > 0.5 )]._xyz );
6054       if ( dist < distMin || projected )
6055       {
6056         _iSeg[ is2nd ] = i;
6057         pProj[ is2nd ] = _offPoints[i]._xyz + ( v01 * uOnSeg ).XYZ();
6058         distMin = dist;
6059       }
6060       else if ( dist > distPrev )
6061       {
6062         if ( ++nbWorse > 3 ) // avoid projection to the middle of a closed EDGE
6063           break;
6064       }
6065       distPrev = dist;
6066       i += di;
6067     }
6068     while ( !projected &&
6069             i >= 0 && i+1 < (int)_offPoints.size() );
6070
6071     if ( !projected )
6072     {
6073       if (( is2nd && _iSeg[1] != _offPoints.size()-2 ) || ( !is2nd && _iSeg[0] != 0 ))
6074       {
6075         _iSeg[0] = 0;
6076         _iSeg[1] = _offPoints.size()-2;
6077         debugMsg( "smoothComplexEdge() failed to project nodes of extreme _LayerEdge's" );
6078         return false;
6079       }
6080     }
6081   }
6082   if ( _iSeg[0] > _iSeg[1] )
6083   {
6084     debugMsg( "smoothComplexEdge() incorrectly projected nodes of extreme _LayerEdge's" );
6085     return false;
6086   }
6087
6088   // adjust length of extreme LE (test viscous_layers_01/B7)
6089   gp_Vec vDiv0( pExtreme[0], pProj[0] );
6090   gp_Vec vDiv1( pExtreme[1], pProj[1] );
6091   double d0 = vDiv0.Magnitude();
6092   double d1 = isProjected[1] ? vDiv1.Magnitude() : 0;
6093   if ( e[0]->Is( _LayerEdge::BLOCKED )) {
6094     if ( e[0]->_normal * vDiv0.XYZ() < 0 ) e[0]->_len += d0;
6095     else                                   e[0]->_len -= d0;
6096   }
6097   if ( e[1]->Is( _LayerEdge::BLOCKED )) {
6098     if ( e[1]->_normal * vDiv1.XYZ() < 0 ) e[1]->_len += d1;
6099     else                                   e[1]->_len -= d1;
6100   }
6101
6102   // ---------------------------------------------------------------------------------
6103   // compute normalized length of the offset segments located between the projections
6104   // ---------------------------------------------------------------------------------
6105
6106   // temporary replace extreme _offPoints by pExtreme
6107   gp_XYZ opXYZ[2] = { _offPoints[ _iSeg[0]   ]._xyz,
6108                       _offPoints[ _iSeg[1]+1 ]._xyz };
6109   _offPoints[ _iSeg[0]   ]._xyz = pExtreme[0].XYZ();
6110   _offPoints[ _iSeg[1]+ 1]._xyz = pExtreme[1].XYZ();
6111
6112   size_t iSeg = 0, nbSeg = _iSeg[1] - _iSeg[0] + 1;
6113   vector< double > len( nbSeg + 1 );
6114   len[ iSeg++ ] = 0;
6115   len[ iSeg++ ] = pProj[ 0 ].Distance( _offPoints[ _iSeg[0]+1 ]._xyz );
6116   for ( size_t i = _iSeg[0]+1; i <= _iSeg[1]; ++i, ++iSeg )
6117   {
6118     len[ iSeg ] = len[ iSeg-1 ] + _offPoints[i].Distance( _offPoints[i+1] );
6119   }
6120   // if ( isProjected[ 1 ])
6121   //   len[ nbSeg ] -= pProj[ 1 ].Distance( _offPoints[ _iSeg[1]+1 ]._xyz );
6122   // else
6123   //   len[ nbSeg ] += pExtreme[ 1 ].Distance( _offPoints[ _iSeg[1]+1 ]._xyz );
6124
6125   double fullLen = len.back() - d0 - d1;
6126   for ( iSeg = 0; iSeg < len.size(); ++iSeg )
6127     len[iSeg] = ( len[iSeg] - d0 ) / fullLen;
6128
6129   // -------------------------------------------------------------
6130   // distribute tgt nodes of _LayerEdge's between the projections
6131   // -------------------------------------------------------------
6132
6133   iSeg = 0;
6134   for ( size_t i = 0; i < _eos.size(); ++i )
6135   {
6136     if ( _eos[i]->Is( _LayerEdge::BLOCKED )) continue;
6137     //if ( !_eos[i]->Is( _LayerEdge::TO_SMOOTH )) continue;
6138     while ( iSeg+2 < len.size() && _leParams[i] > len[ iSeg+1 ] )
6139       iSeg++;
6140     double r = ( _leParams[i] - len[ iSeg ]) / ( len[ iSeg+1 ] - len[ iSeg ]);
6141     gp_XYZ p = ( _offPoints[ iSeg + _iSeg[0]     ]._xyz * ( 1 - r ) +
6142                  _offPoints[ iSeg + _iSeg[0] + 1 ]._xyz * r );
6143
6144     if ( surface.IsNull() )
6145     {
6146       _eos[i]->_pos.back() = p;
6147     }
6148     else // project a new node position to a FACE
6149     {
6150       gp_Pnt2d uv ( _eos[i]->_pos.back().X(), _eos[i]->_pos.back().Y() );
6151       gp_Pnt2d uv2( surface->NextValueOfUV( uv, p, fTol ));
6152
6153       p = surface->Value( uv2 ).XYZ();
6154       _eos[i]->_pos.back().SetCoord( uv2.X(), uv2.Y(), 0 );
6155     }
6156     SMDS_MeshNode* tgtNode = const_cast<SMDS_MeshNode*>( _eos[i]->_nodes.back() );
6157     tgtNode->setXYZ( p.X(), p.Y(), p.Z() );
6158     dumpMove( tgtNode );
6159   }
6160
6161   _offPoints[ _iSeg[0]   ]._xyz = opXYZ[0];
6162   _offPoints[ _iSeg[1]+1 ]._xyz = opXYZ[1];
6163
6164   return true;
6165 }
6166
6167 //================================================================================
6168 /*!
6169  * \brief Prepare for smoothing
6170  */
6171 //================================================================================
6172
6173 void _Smoother1D::prepare(_SolidData& data)
6174 {
6175   const TopoDS_Edge& E = TopoDS::Edge( _eos._shape );
6176   _curveLen = SMESH_Algo::EdgeLength( E );
6177
6178   // sort _LayerEdge's by position on the EDGE
6179   data.SortOnEdge( E, _eos._edges );
6180
6181   // compute normalized param of _eos._edges on EDGE
6182   _leParams.resize( _eos._edges.size() + 1 );
6183   {
6184     double curLen;
6185     gp_Pnt pPrev = SMESH_TNodeXYZ( getLEdgeOnV( 0 )->_nodes[0] );
6186     _leParams[0] = 0;
6187     for ( size_t i = 0; i < _eos._edges.size(); ++i )
6188     {
6189       gp_Pnt p       = SMESH_TNodeXYZ( _eos._edges[i]->_nodes[0] );
6190       curLen         = p.Distance( pPrev );
6191       _leParams[i+1] = _leParams[i] + curLen;
6192       pPrev          = p;
6193     }
6194     double fullLen = _leParams.back() + pPrev.Distance( SMESH_TNodeXYZ( getLEdgeOnV(1)->_nodes[0]));
6195     for ( size_t i = 0; i < _leParams.size()-1; ++i )
6196       _leParams[i] = _leParams[i+1] / fullLen;
6197     _leParams.back() = 1.;
6198   }
6199
6200   _LayerEdge* leOnV[2] = { getLEdgeOnV(0), getLEdgeOnV(1) };
6201
6202   // get cosin to use in findEdgesToSmooth()
6203   _edgeDir[0] = getEdgeDir( E, leOnV[0]->_nodes[0], data.GetHelper() );
6204   _edgeDir[1] = getEdgeDir( E, leOnV[1]->_nodes[0], data.GetHelper() );
6205   _leOnV[0]._cosin = Abs( leOnV[0]->_cosin );
6206   _leOnV[1]._cosin = Abs( leOnV[1]->_cosin );
6207   if ( _eos._sWOL.IsNull() ) // 3D
6208     for ( int iEnd = 0; iEnd < 2; ++iEnd )
6209       _leOnV[iEnd]._cosin = Abs( _edgeDir[iEnd].Normalized() * leOnV[iEnd]->_normal );
6210
6211   if ( isAnalytic() )
6212     return;
6213
6214   // divide E to have offset segments with low deflection
6215   BRepAdaptor_Curve c3dAdaptor( E );
6216   const double curDeflect = 0.1; //0.01; // Curvature deflection == |p1p2]*sin(p1p2,p1pM)
6217   const double angDeflect = 0.1; //0.09; // Angular deflection == sin(p1pM,pMp2)
6218   GCPnts_TangentialDeflection discret(c3dAdaptor, angDeflect, curDeflect);
6219   if ( discret.NbPoints() <= 2 )
6220   {
6221     _anaCurve = new Geom_Line( gp::OX() ); // only type does matter
6222     return;
6223   }
6224
6225   const double u0 = c3dAdaptor.FirstParameter();
6226   gp_Pnt p; gp_Vec tangent;
6227   if ( discret.NbPoints() >= (int) _eos.size() + 2 )
6228   {
6229     _offPoints.resize( discret.NbPoints() );
6230     for ( size_t i = 0; i < _offPoints.size(); i++ )
6231     {
6232       double u = discret.Parameter( i+1 );
6233       c3dAdaptor.D1( u, p, tangent );
6234       _offPoints[i]._xyz     = p.XYZ();
6235       _offPoints[i]._edgeDir = tangent.XYZ();
6236       _offPoints[i]._param = GCPnts_AbscissaPoint::Length( c3dAdaptor, u0, u ) / _curveLen;
6237     }
6238   }
6239   else
6240   {
6241     std::vector< double > params( _eos.size() + 2 );
6242
6243     params[0]     = data.GetHelper().GetNodeU( E, leOnV[0]->_nodes[0] );
6244     params.back() = data.GetHelper().GetNodeU( E, leOnV[1]->_nodes[0] );
6245     for ( size_t i = 0; i < _eos.size(); i++ )
6246       params[i+1] = data.GetHelper().GetNodeU( E, _eos[i]->_nodes[0] );
6247
6248     if ( params[1] > params[ _eos.size() ] )
6249       std::reverse( params.begin() + 1, params.end() - 1 );
6250
6251     _offPoints.resize( _eos.size() + 2 );
6252     for ( size_t i = 0; i < _offPoints.size(); i++ )
6253     {
6254       const double u = params[i];
6255       c3dAdaptor.D1( u, p, tangent );
6256       _offPoints[i]._xyz     = p.XYZ();
6257       _offPoints[i]._edgeDir = tangent.XYZ();
6258       _offPoints[i]._param = GCPnts_AbscissaPoint::Length( c3dAdaptor, u0, u ) / _curveLen;
6259     }
6260   }
6261
6262   // set _2edges
6263   _offPoints    [0]._2edges.set( &_leOnV[0], &_leOnV[0], 0.5, 0.5 );
6264   _offPoints.back()._2edges.set( &_leOnV[1], &_leOnV[1], 0.5, 0.5 );
6265   _2NearEdges tmp2edges;
6266   tmp2edges._edges[1] = _eos._edges[0];
6267   _leOnV[0]._2neibors = & tmp2edges;
6268   _leOnV[0]._nodes    = leOnV[0]->_nodes;
6269   _leOnV[1]._nodes    = leOnV[1]->_nodes;
6270   _LayerEdge* eNext, *ePrev = & _leOnV[0];
6271   for ( size_t iLE = 0, i = 1; i < _offPoints.size()-1; i++ )
6272   {
6273     // find _LayerEdge's located before and after an offset point
6274     // (_eos._edges[ iLE ] is next after ePrev)
6275     while ( iLE < _eos._edges.size() && _offPoints[i]._param > _leParams[ iLE ] )
6276       ePrev = _eos._edges[ iLE++ ];
6277     eNext = ePrev->_2neibors->_edges[1];
6278
6279     gp_Pnt p0 = SMESH_TNodeXYZ( ePrev->_nodes[0] );
6280     gp_Pnt p1 = SMESH_TNodeXYZ( eNext->_nodes[0] );
6281     double  r = p0.Distance( _offPoints[i]._xyz ) / p0.Distance( p1 );
6282     _offPoints[i]._2edges.set( ePrev, eNext, 1-r, r );
6283   }
6284
6285   // replace _LayerEdge's on VERTEX by _leOnV in _offPoints._2edges
6286   for ( size_t i = 0; i < _offPoints.size(); i++ )
6287     if ( _offPoints[i]._2edges._edges[0] == leOnV[0] )
6288       _offPoints[i]._2edges._edges[0] = & _leOnV[0];
6289     else break;
6290   for ( size_t i = _offPoints.size()-1; i > 0; i-- )
6291     if ( _offPoints[i]._2edges._edges[1] == leOnV[1] )
6292       _offPoints[i]._2edges._edges[1] = & _leOnV[1];
6293     else break;
6294
6295   // set _normal of _leOnV[0] and _leOnV[1] to be normal to the EDGE
6296
6297   int iLBO = _offPoints.size() - 2; // last but one
6298
6299   if ( leOnV[ 0 ]->Is( _LayerEdge::MULTI_NORMAL ))
6300     _leOnV[ 0 ]._normal = getNormalNormal( _eos._edges[1]->_normal, _edgeDir[0] );
6301   else
6302     _leOnV[ 0 ]._normal = getNormalNormal( leOnV[0]->_normal,       _edgeDir[0] );
6303   if ( leOnV[ 1 ]->Is( _LayerEdge::MULTI_NORMAL ))
6304     _leOnV[ 1 ]._normal = getNormalNormal( _eos._edges.back()->_normal, _edgeDir[1] );
6305   else
6306     _leOnV[ 1 ]._normal = getNormalNormal( leOnV[1]->_normal,           _edgeDir[1] );
6307   _leOnV[ 0 ]._len = 0;
6308   _leOnV[ 1 ]._len = 0;
6309   _leOnV[ 0 ]._lenFactor = _offPoints[1   ]._2edges._edges[1]->_lenFactor;
6310   _leOnV[ 1 ]._lenFactor = _offPoints[iLBO]._2edges._edges[0]->_lenFactor;
6311
6312   _iSeg[0] = 0;
6313   _iSeg[1] = _offPoints.size()-2;
6314
6315   // initialize OffPnt::_len
6316   for ( size_t i = 0; i < _offPoints.size(); ++i )
6317     _offPoints[i]._len = 0;
6318
6319   if ( _eos._edges[0]->NbSteps() > 1 ) // already inflated several times, init _xyz
6320   {
6321     _leOnV[0]._len = leOnV[0]->_len;
6322     _leOnV[1]._len = leOnV[1]->_len;
6323     for ( size_t i = 0; i < _offPoints.size(); i++ )
6324     {
6325       _LayerEdge*  e0 = _offPoints[i]._2edges._edges[0];
6326       _LayerEdge*  e1 = _offPoints[i]._2edges._edges[1];
6327       const double w0 = _offPoints[i]._2edges._wgt[0];
6328       const double w1 = _offPoints[i]._2edges._wgt[1];
6329       double  avgLen  = ( e0->_len * w0 + e1->_len * w1 );
6330       gp_XYZ  avgXYZ  = ( SMESH_TNodeXYZ( e0->_nodes.back() ) * w0 +
6331                           SMESH_TNodeXYZ( e1->_nodes.back() ) * w1 );
6332       _offPoints[i]._xyz = avgXYZ;
6333       _offPoints[i]._len = avgLen;
6334     }
6335   }
6336 }
6337
6338 //================================================================================
6339 /*!
6340  * \brief return _normal of _leOnV[is2nd] normal to the EDGE
6341  */
6342 //================================================================================
6343
6344 gp_XYZ _Smoother1D::getNormalNormal( const gp_XYZ & normal,
6345                                      const gp_XYZ&  edgeDir)
6346 {
6347   gp_XYZ cross = normal ^ edgeDir;
6348   gp_XYZ  norm = edgeDir ^ cross;
6349   double  size = norm.Modulus();
6350
6351   // if ( size == 0 ) // MULTI_NORMAL _LayerEdge
6352   //   return gp_XYZ( 1e-100, 1e-100, 1e-100 );
6353
6354   return norm / size;
6355 }
6356
6357 //================================================================================
6358 /*!
6359  * \brief Writes a script creating a mesh composed of _offPoints
6360  */
6361 //================================================================================
6362
6363 void _Smoother1D::offPointsToPython() const
6364 {
6365   const char* fname = "/tmp/offPoints.py";
6366   cout << "execfile('"<<fname<<"')"<<endl;
6367   ofstream py(fname);
6368   py << "import SMESH" << endl
6369      << "from salome.smesh import smeshBuilder" << endl
6370      << "smesh  = smeshBuilder.New(salome.myStudy)" << endl
6371      << "mesh   = smesh.Mesh( 'offPoints' )"<<endl;
6372   for ( size_t i = 0; i < _offPoints.size(); i++ )
6373   {
6374     py << "mesh.AddNode( "
6375        << _offPoints[i]._xyz.X() << ", "
6376        << _offPoints[i]._xyz.Y() << ", "
6377        << _offPoints[i]._xyz.Z() << " )" << endl;
6378   }
6379 }
6380
6381 //================================================================================
6382 /*!
6383  * \brief Sort _LayerEdge's by a parameter on a given EDGE
6384  */
6385 //================================================================================
6386
6387 void _SolidData::SortOnEdge( const TopoDS_Edge&     E,
6388                              vector< _LayerEdge* >& edges)
6389 {
6390   map< double, _LayerEdge* > u2edge;
6391   for ( size_t i = 0; i < edges.size(); ++i )
6392     u2edge.insert( u2edge.end(),
6393                    make_pair( _helper->GetNodeU( E, edges[i]->_nodes[0] ), edges[i] ));
6394
6395   ASSERT( u2edge.size() == edges.size() );
6396   map< double, _LayerEdge* >::iterator u2e = u2edge.begin();
6397   for ( size_t i = 0; i < edges.size(); ++i, ++u2e )
6398     edges[i] = u2e->second;
6399
6400   Sort2NeiborsOnEdge( edges );
6401 }
6402
6403 //================================================================================
6404 /*!
6405  * \brief Set _2neibors according to the order of _LayerEdge on EDGE
6406  */
6407 //================================================================================
6408
6409 void _SolidData::Sort2NeiborsOnEdge( vector< _LayerEdge* >& edges )
6410 {
6411   if ( edges.size() < 2 || !edges[0]->_2neibors ) return;
6412
6413   for ( size_t i = 0; i < edges.size()-1; ++i )
6414     if ( edges[i]->_2neibors->tgtNode(1) != edges[i+1]->_nodes.back() )
6415       edges[i]->_2neibors->reverse();
6416
6417   const size_t iLast = edges.size() - 1;
6418   if ( edges.size() > 1 &&
6419        edges[iLast]->_2neibors->tgtNode(0) != edges[iLast-1]->_nodes.back() )
6420     edges[iLast]->_2neibors->reverse();
6421 }
6422
6423 //================================================================================
6424 /*!
6425  * \brief Return _EdgesOnShape* corresponding to the shape
6426  */
6427 //================================================================================
6428
6429 _EdgesOnShape* _SolidData::GetShapeEdges(const TGeomID shapeID )
6430 {
6431   if ( shapeID < (int)_edgesOnShape.size() &&
6432        _edgesOnShape[ shapeID ]._shapeID == shapeID )
6433     return _edgesOnShape[ shapeID ]._subMesh ? & _edgesOnShape[ shapeID ] : 0;
6434
6435   for ( size_t i = 0; i < _edgesOnShape.size(); ++i )
6436     if ( _edgesOnShape[i]._shapeID == shapeID )
6437       return _edgesOnShape[i]._subMesh ? & _edgesOnShape[i] : 0;
6438
6439   return 0;
6440 }
6441
6442 //================================================================================
6443 /*!
6444  * \brief Return _EdgesOnShape* corresponding to the shape
6445  */
6446 //================================================================================
6447
6448 _EdgesOnShape* _SolidData::GetShapeEdges(const TopoDS_Shape& shape )
6449 {
6450   SMESHDS_Mesh* meshDS = _proxyMesh->GetMesh()->GetMeshDS();
6451   return GetShapeEdges( meshDS->ShapeToIndex( shape ));
6452 }
6453
6454 //================================================================================
6455 /*!
6456  * \brief Prepare data of the _LayerEdge for smoothing on FACE
6457  */
6458 //================================================================================
6459
6460 void _SolidData::PrepareEdgesToSmoothOnFace( _EdgesOnShape* eos, bool substituteSrcNodes )
6461 {
6462   SMESH_MesherHelper helper( *_proxyMesh->GetMesh() );
6463
6464   set< TGeomID > vertices;
6465   TopoDS_Face F;
6466   if ( eos->ShapeType() == TopAbs_FACE )
6467   {
6468     // check FACE concavity and get concave VERTEXes
6469     F = TopoDS::Face( eos->_shape );
6470     if ( isConcave( F, helper, &vertices ))
6471       _concaveFaces.insert( eos->_shapeID );
6472
6473     // set eos._eosConcaVer
6474     eos->_eosConcaVer.clear();
6475     eos->_eosConcaVer.reserve( vertices.size() );
6476     for ( set< TGeomID >::iterator v = vertices.begin(); v != vertices.end(); ++v )
6477     {
6478       _EdgesOnShape* eov = GetShapeEdges( *v );
6479       if ( eov && eov->_edges.size() == 1 )
6480       {
6481         eos->_eosConcaVer.push_back( eov );
6482         for ( size_t i = 0; i < eov->_edges[0]->_neibors.size(); ++i )
6483           eov->_edges[0]->_neibors[i]->Set( _LayerEdge::DIFFICULT );
6484       }
6485     }
6486
6487     // SetSmooLen() to _LayerEdge's on FACE
6488     // for ( size_t i = 0; i < eos->_edges.size(); ++i )
6489     // {
6490     //   eos->_edges[i]->SetSmooLen( Precision::Infinite() );
6491     // }
6492     // SMESH_subMeshIteratorPtr smIt = eos->_subMesh->getDependsOnIterator(/*includeSelf=*/false);
6493     // while ( smIt->more() ) // loop on sub-shapes of the FACE
6494     // {
6495     //   _EdgesOnShape* eoe = GetShapeEdges( smIt->next()->GetId() );
6496     //   if ( !eoe ) continue;
6497
6498     //   vector<_LayerEdge*>& eE = eoe->_edges;
6499     //   for ( size_t iE = 0; iE < eE.size(); ++iE ) // loop on _LayerEdge's on EDGE or VERTEX
6500     //   {
6501     //     if ( eE[iE]->_cosin <= theMinSmoothCosin )
6502     //       continue;
6503
6504     //     SMDS_ElemIteratorPtr segIt = eE[iE]->_nodes[0]->GetInverseElementIterator(SMDSAbs_Edge);
6505     //     while ( segIt->more() )
6506     //     {
6507     //       const SMDS_MeshElement* seg = segIt->next();
6508     //       if ( !eos->_subMesh->DependsOn( seg->getshapeId() ))
6509     //         continue;
6510     //       if ( seg->GetNode(0) != eE[iE]->_nodes[0] )
6511     //         continue; // not to check a seg twice
6512     //       for ( size_t iN = 0; iN < eE[iE]->_neibors.size(); ++iN )
6513     //       {
6514     //         _LayerEdge* eN = eE[iE]->_neibors[iN];
6515     //         if ( eN->_nodes[0]->getshapeId() != eos->_shapeID )
6516     //           continue;
6517     //         double dist    = SMESH_MeshAlgos::GetDistance( seg, SMESH_TNodeXYZ( eN->_nodes[0] ));
6518     //         double smooLen = getSmoothingThickness( eE[iE]->_cosin, dist );
6519     //         eN->SetSmooLen( Min( smooLen, eN->GetSmooLen() ));
6520     //         eN->Set( _LayerEdge::NEAR_BOUNDARY );
6521     //       }
6522     //     }
6523     //   }
6524     // }
6525   } // if ( eos->ShapeType() == TopAbs_FACE )
6526
6527   for ( size_t i = 0; i < eos->_edges.size(); ++i )
6528   {
6529     eos->_edges[i]->_smooFunction = 0;
6530     eos->_edges[i]->Set( _LayerEdge::TO_SMOOTH );
6531   }
6532   bool isCurved = false;
6533   for ( size_t i = 0; i < eos->_edges.size(); ++i )
6534   {
6535     _LayerEdge* edge = eos->_edges[i];
6536
6537     // get simplices sorted
6538     _Simplex::SortSimplices( edge->_simplices );
6539
6540     // smoothing function
6541     edge->ChooseSmooFunction( vertices, _n2eMap );
6542
6543     // set _curvature
6544     double avgNormProj = 0, avgLen = 0;
6545     for ( size_t iS = 0; iS < edge->_simplices.size(); ++iS )
6546     {
6547       _Simplex& s = edge->_simplices[iS];
6548
6549       gp_XYZ  vec = edge->_pos.back() - SMESH_TNodeXYZ( s._nPrev );
6550       avgNormProj += edge->_normal * vec;
6551       avgLen      += vec.Modulus();
6552       if ( substituteSrcNodes )
6553       {
6554         s._nNext = _n2eMap[ s._nNext ]->_nodes.back();
6555         s._nPrev = _n2eMap[ s._nPrev ]->_nodes.back();
6556       }
6557     }
6558     avgNormProj /= edge->_simplices.size();
6559     avgLen      /= edge->_simplices.size();
6560     if (( edge->_curvature = _Curvature::New( avgNormProj, avgLen )))
6561     {
6562       edge->Set( _LayerEdge::SMOOTHED_C1 );
6563       isCurved = true;
6564       SMDS_FacePositionPtr fPos = edge->_nodes[0]->GetPosition();
6565       if ( !fPos )
6566         for ( size_t iS = 0; iS < edge->_simplices.size()  &&  !fPos; ++iS )
6567           fPos = edge->_simplices[iS]._nPrev->GetPosition();
6568       if ( fPos )
6569         edge->_curvature->_uv.SetCoord( fPos->GetUParameter(), fPos->GetVParameter() );
6570     }
6571   }
6572
6573   // prepare for putOnOffsetSurface()
6574   if (( eos->ShapeType() == TopAbs_FACE ) &&
6575       ( isCurved || !eos->_eosConcaVer.empty() ))
6576   {
6577     eos->_offsetSurf = helper.GetSurface( TopoDS::Face( eos->_shape ));
6578     eos->_edgeForOffset = 0;
6579
6580     double maxCosin = -1;
6581     for ( TopExp_Explorer eExp( eos->_shape, TopAbs_EDGE ); eExp.More(); eExp.Next() )
6582     {
6583       _EdgesOnShape* eoe = GetShapeEdges( eExp.Current() );
6584       if ( !eoe || eoe->_edges.empty() ) continue;
6585
6586       vector<_LayerEdge*>& eE = eoe->_edges;
6587       _LayerEdge* e = eE[ eE.size() / 2 ];
6588       if ( e->_cosin > maxCosin )
6589       {
6590         eos->_edgeForOffset = e;
6591         maxCosin = e->_cosin;
6592       }
6593     }
6594   }
6595 }
6596
6597 //================================================================================
6598 /*!
6599  * \brief Add faces for smoothing
6600  */
6601 //================================================================================
6602
6603 void _SolidData::AddShapesToSmooth( const set< _EdgesOnShape* >& eosToSmooth,
6604                                     const set< _EdgesOnShape* >* edgesNoAnaSmooth )
6605 {
6606   set< _EdgesOnShape * >::const_iterator eos = eosToSmooth.begin();
6607   for ( ; eos != eosToSmooth.end(); ++eos )
6608   {
6609     if ( !*eos || (*eos)->_toSmooth ) continue;
6610
6611     (*eos)->_toSmooth = true;
6612
6613     if ( (*eos)->ShapeType() == TopAbs_FACE )
6614     {
6615       PrepareEdgesToSmoothOnFace( *eos, /*substituteSrcNodes=*/false );
6616       (*eos)->_toSmooth = true;
6617     }
6618   }
6619
6620   // avoid _Smoother1D::smoothAnalyticEdge() of edgesNoAnaSmooth
6621   if ( edgesNoAnaSmooth )
6622     for ( eos = edgesNoAnaSmooth->begin(); eos != edgesNoAnaSmooth->end(); ++eos )
6623     {
6624       if ( (*eos)->_edgeSmoother )
6625         (*eos)->_edgeSmoother->_anaCurve.Nullify();
6626     }
6627 }
6628
6629 //================================================================================
6630 /*!
6631  * \brief Limit _LayerEdge::_maxLen according to local curvature
6632  */
6633 //================================================================================
6634
6635 void _ViscousBuilder::limitMaxLenByCurvature( _SolidData& data, SMESH_MesherHelper& helper )
6636 {
6637   // find intersection of neighbor _LayerEdge's to limit _maxLen
6638   // according to local curvature (IPAL52648)
6639
6640   // This method must be called after findCollisionEdges() where _LayerEdge's
6641   // get _lenFactor initialized in the case of eos._hyp.IsOffsetMethod()
6642
6643   for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
6644   {
6645     _EdgesOnShape& eosI = data._edgesOnShape[iS];
6646     if ( eosI._edges.empty() ) continue;
6647     if ( !eosI._hyp.ToSmooth() )
6648     {
6649       for ( size_t i = 0; i < eosI._edges.size(); ++i )
6650       {
6651         _LayerEdge* eI = eosI._edges[i];
6652         for ( size_t iN = 0; iN < eI->_neibors.size(); ++iN )
6653         {
6654           _LayerEdge* eN = eI->_neibors[iN];
6655           if ( eI->_nodes[0]->GetID() < eN->_nodes[0]->GetID() ) // treat this pair once
6656           {
6657             _EdgesOnShape* eosN = data.GetShapeEdges( eN );
6658             limitMaxLenByCurvature( eI, eN, eosI, *eosN, eosI._hyp.ToSmooth() );
6659           }
6660         }
6661       }
6662     }
6663     else if ( eosI.ShapeType() == TopAbs_EDGE )
6664     {
6665       const TopoDS_Edge& E = TopoDS::Edge( eosI._shape );
6666       if ( SMESH_Algo::IsStraight( E, /*degenResult=*/true )) continue;
6667
6668       _LayerEdge* e0 = eosI._edges[0];
6669       for ( size_t i = 1; i < eosI._edges.size(); ++i )
6670       {
6671         _LayerEdge* eI = eosI._edges[i];
6672         limitMaxLenByCurvature( eI, e0, eosI, eosI, eosI._hyp.ToSmooth() );
6673         e0 = eI;
6674       }
6675     }
6676   }
6677 }
6678
6679 //================================================================================
6680 /*!
6681  * \brief Limit _LayerEdge::_maxLen according to local curvature
6682  */
6683 //================================================================================
6684
6685 void _ViscousBuilder::limitMaxLenByCurvature( _LayerEdge*    e1,
6686                                               _LayerEdge*    e2,
6687                                               _EdgesOnShape& eos1,
6688                                               _EdgesOnShape& eos2,
6689                                               const bool     isSmoothable )
6690 {
6691   if (( e1->_nodes[0]->GetPosition()->GetDim() !=
6692         e2->_nodes[0]->GetPosition()->GetDim() ) &&
6693       ( e1->_cosin < 0.75 ))
6694     return; // angle > 90 deg at e1
6695
6696   gp_XYZ plnNorm = e1->_normal ^ e2->_normal;
6697   double norSize = plnNorm.SquareModulus();
6698   if ( norSize < std::numeric_limits<double>::min() )
6699     return; // parallel normals
6700
6701   // find closest points of skew _LayerEdge's
6702   SMESH_TNodeXYZ src1( e1->_nodes[0] ), src2( e2->_nodes[0] );
6703   gp_XYZ dir12 = src2 - src1;
6704   gp_XYZ perp1 = e1->_normal ^ plnNorm;
6705   gp_XYZ perp2 = e2->_normal ^ plnNorm;
6706   double  dot1 = perp2 * e1->_normal;
6707   double  dot2 = perp1 * e2->_normal;
6708   double    u1 =   ( perp2 * dir12 ) / dot1;
6709   double    u2 = - ( perp1 * dir12 ) / dot2;
6710   if ( u1 > 0 && u2 > 0 )
6711   {
6712     double ovl = ( u1 * e1->_normal * dir12 -
6713                    u2 * e2->_normal * dir12 ) / dir12.SquareModulus();
6714     if ( ovl > theSmoothThickToElemSizeRatio )
6715     {
6716       const double coef = 0.75;
6717       e1->SetMaxLen( Min( e1->_maxLen, coef * u1 / e1->_lenFactor ));
6718       e2->SetMaxLen( Min( e2->_maxLen, coef * u2 / e2->_lenFactor ));
6719     }
6720   }
6721 }
6722
6723 //================================================================================
6724 /*!
6725  * \brief Fill data._collisionEdges
6726  */
6727 //================================================================================
6728
6729 void _ViscousBuilder::findCollisionEdges( _SolidData& data, SMESH_MesherHelper& helper )
6730 {
6731   data._collisionEdges.clear();
6732
6733   // set the full thickness of the layers to LEs
6734   for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
6735   {
6736     _EdgesOnShape& eos = data._edgesOnShape[iS];
6737     if ( eos._edges.empty() ) continue;
6738     if ( eos.ShapeType() != TopAbs_EDGE && eos.ShapeType() != TopAbs_VERTEX ) continue;
6739     if ( !eos._sWOL.IsNull() ) continue; // PAL23566
6740
6741     for ( size_t i = 0; i < eos._edges.size(); ++i )
6742     {
6743       if ( eos._edges[i]->Is( _LayerEdge::BLOCKED )) continue;
6744       double maxLen = eos._edges[i]->_maxLen;
6745       eos._edges[i]->_maxLen = Precision::Infinite(); // avoid blocking
6746       eos._edges[i]->SetNewLength( 1.5 * maxLen, eos, helper );
6747       eos._edges[i]->_maxLen = maxLen;
6748     }
6749   }
6750
6751   // make temporary quadrangles got by extrusion of
6752   // mesh edges along _LayerEdge._normal's
6753
6754   vector< const SMDS_MeshElement* > tmpFaces;
6755
6756   for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
6757   {
6758     _EdgesOnShape& eos = data._edgesOnShape[ iS ];
6759     if ( eos.ShapeType() != TopAbs_EDGE )
6760       continue;
6761     if ( eos._edges.empty() )
6762     {
6763       _LayerEdge* edge[2] = { 0, 0 }; // LE of 2 VERTEX'es
6764       SMESH_subMeshIteratorPtr smIt = eos._subMesh->getDependsOnIterator(/*includeSelf=*/false);
6765       while ( smIt->more() )
6766         if ( _EdgesOnShape* eov = data.GetShapeEdges( smIt->next()->GetId() ))
6767           if ( eov->_edges.size() == 1 )
6768             edge[ bool( edge[0]) ] = eov->_edges[0];
6769
6770       if ( edge[1] )
6771       {
6772         _TmpMeshFaceOnEdge* f = new _TmpMeshFaceOnEdge( edge[0], edge[1], --_tmpFaceID );
6773         tmpFaces.push_back( f );
6774       }
6775     }
6776     for ( size_t i = 0; i < eos._edges.size(); ++i )
6777     {
6778       _LayerEdge* edge = eos._edges[i];
6779       for ( int j = 0; j < 2; ++j ) // loop on _2NearEdges
6780       {
6781         const SMDS_MeshNode* src2 = edge->_2neibors->srcNode(j);
6782         if ( src2->GetPosition()->GetDim() > 0 &&
6783              src2->GetID() < edge->_nodes[0]->GetID() )
6784           continue; // avoid using same segment twice
6785
6786         // a _LayerEdge containing tgt2
6787         _LayerEdge* neiborEdge = edge->_2neibors->_edges[j];
6788
6789         _TmpMeshFaceOnEdge* f = new _TmpMeshFaceOnEdge( edge, neiborEdge, --_tmpFaceID );
6790         tmpFaces.push_back( f );
6791       }
6792     }
6793   }
6794
6795   // Find _LayerEdge's intersecting tmpFaces.
6796
6797   SMDS_ElemIteratorPtr fIt( new SMDS_ElementVectorIterator( tmpFaces.begin(),
6798                                                             tmpFaces.end()));
6799   SMESHUtils::Deleter<SMESH_ElementSearcher> searcher
6800     ( SMESH_MeshAlgos::GetElementSearcher( *getMeshDS(), fIt ));
6801
6802   double dist1, dist2, segLen, eps = 0.5;
6803   _CollisionEdges collEdges;
6804   vector< const SMDS_MeshElement* > suspectFaces;
6805   const double angle45 = Cos( 45. * M_PI / 180. );
6806
6807   for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
6808   {
6809     _EdgesOnShape& eos = data._edgesOnShape[ iS ];
6810     if ( eos.ShapeType() == TopAbs_FACE || !eos._sWOL.IsNull() )
6811       continue;
6812     // find sub-shapes whose VL can influence VL on eos
6813     set< TGeomID > neighborShapes;
6814     PShapeIteratorPtr fIt = helper.GetAncestors( eos._shape, *_mesh, TopAbs_FACE );
6815     while ( const TopoDS_Shape* face = fIt->next() )
6816     {
6817       TGeomID faceID = getMeshDS()->ShapeToIndex( *face );
6818       if ( _EdgesOnShape* eof = data.GetShapeEdges( faceID ))
6819       {
6820         SMESH_subMeshIteratorPtr subIt = eof->_subMesh->getDependsOnIterator(/*includeSelf=*/false);
6821         while ( subIt->more() )
6822           neighborShapes.insert( subIt->next()->GetId() );
6823       }
6824     }
6825     if ( eos.ShapeType() == TopAbs_VERTEX )
6826     {
6827       PShapeIteratorPtr eIt = helper.GetAncestors( eos._shape, *_mesh, TopAbs_EDGE );
6828       while ( const TopoDS_Shape* edge = eIt->next() )
6829         neighborShapes.erase( getMeshDS()->ShapeToIndex( *edge ));
6830     }
6831     // find intersecting _LayerEdge's
6832     for ( size_t i = 0; i < eos._edges.size(); ++i )
6833     {
6834       if ( eos._edges[i]->Is( _LayerEdge::MULTI_NORMAL )) continue;
6835       _LayerEdge*   edge = eos._edges[i];
6836       gp_Ax1 lastSegment = edge->LastSegment( segLen, eos );
6837       segLen *= 1.2;
6838
6839       gp_Vec eSegDir0, eSegDir1;
6840       if ( edge->IsOnEdge() )
6841       {
6842         SMESH_TNodeXYZ eP( edge->_nodes[0] );
6843         eSegDir0 = SMESH_TNodeXYZ( edge->_2neibors->srcNode(0) ) - eP;
6844         eSegDir1 = SMESH_TNodeXYZ( edge->_2neibors->srcNode(1) ) - eP;
6845       }
6846       suspectFaces.clear();
6847       searcher->GetElementsInSphere( SMESH_TNodeXYZ( edge->_nodes.back()), edge->_len * 2,
6848                                      SMDSAbs_Face, suspectFaces );
6849       collEdges._intEdges.clear();
6850       for ( size_t j = 0 ; j < suspectFaces.size(); ++j )
6851       {
6852         const _TmpMeshFaceOnEdge* f = (const _TmpMeshFaceOnEdge*) suspectFaces[j];
6853         if ( f->_le1 == edge || f->_le2 == edge ) continue;
6854         if ( !neighborShapes.count( f->_le1->_nodes[0]->getshapeId() )) continue;
6855         if ( !neighborShapes.count( f->_le2->_nodes[0]->getshapeId() )) continue;
6856         if ( edge->IsOnEdge() ) {
6857           if ( edge->_2neibors->include( f->_le1 ) ||
6858                edge->_2neibors->include( f->_le2 )) continue;
6859         }
6860         else {
6861           if (( f->_le1->IsOnEdge() && f->_le1->_2neibors->include( edge )) ||
6862               ( f->_le2->IsOnEdge() && f->_le2->_2neibors->include( edge )))  continue;
6863         }
6864         dist1 = dist2 = Precision::Infinite();
6865         if ( !edge->SegTriaInter( lastSegment, f->n(0), f->n(1), f->n(2), dist1, eps ))
6866           dist1 = Precision::Infinite();
6867         if ( !edge->SegTriaInter( lastSegment, f->n(3), f->n(2), f->n(0), dist2, eps ))
6868           dist2 = Precision::Infinite();
6869         if (( dist1 > segLen ) && ( dist2 > segLen ))
6870           continue;
6871
6872         if ( edge->IsOnEdge() )
6873         {
6874           // skip perpendicular EDGEs
6875           gp_Vec fSegDir  = SMESH_TNodeXYZ( f->n(0) ) - SMESH_TNodeXYZ( f->n(3) );
6876           bool isParallel = ( isLessAngle( eSegDir0, fSegDir, angle45 ) ||
6877                               isLessAngle( eSegDir1, fSegDir, angle45 ) ||
6878                               isLessAngle( eSegDir0, fSegDir.Reversed(), angle45 ) ||
6879                               isLessAngle( eSegDir1, fSegDir.Reversed(), angle45 ));
6880           if ( !isParallel )
6881             continue;
6882         }
6883
6884         // either limit inflation of edges or remember them for updating _normal
6885         // double dot = edge->_normal * f->GetDir();
6886         // if ( dot > 0.1 )
6887         {
6888           collEdges._intEdges.push_back( f->_le1 );
6889           collEdges._intEdges.push_back( f->_le2 );
6890         }
6891         // else
6892         // {
6893         //   double shortLen = 0.75 * ( Min( dist1, dist2 ) / edge->_lenFactor );
6894         //   edge->SetMaxLen( Min( shortLen, edge->_maxLen ));
6895         // }
6896       }
6897
6898       if ( !collEdges._intEdges.empty() )
6899       {
6900         collEdges._edge = edge;
6901         data._collisionEdges.push_back( collEdges );
6902       }
6903     }
6904   }
6905
6906   for ( size_t i = 0 ; i < tmpFaces.size(); ++i )
6907     delete tmpFaces[i];
6908
6909   // restore the zero thickness
6910   for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
6911   {
6912     _EdgesOnShape& eos = data._edgesOnShape[iS];
6913     if ( eos._edges.empty() ) continue;
6914     if ( eos.ShapeType() != TopAbs_EDGE && eos.ShapeType() != TopAbs_VERTEX ) continue;
6915
6916     for ( size_t i = 0; i < eos._edges.size(); ++i )
6917     {
6918       eos._edges[i]->InvalidateStep( 1, eos );
6919       eos._edges[i]->_len = 0;
6920     }
6921   }
6922 }
6923
6924 //================================================================================
6925 /*!
6926  * \brief Find _LayerEdge's located on boundary of a convex FACE whose normal
6927  *        will be updated at each inflation step
6928  */
6929 //================================================================================
6930
6931 void _ViscousBuilder::findEdgesToUpdateNormalNearConvexFace( _ConvexFace &       convFace,
6932                                                              _SolidData&         data,
6933                                                              SMESH_MesherHelper& helper )
6934 {
6935   const TGeomID convFaceID = getMeshDS()->ShapeToIndex( convFace._face );
6936   const double       preci = BRep_Tool::Tolerance( convFace._face );
6937   Handle(ShapeAnalysis_Surface) surface = helper.GetSurface( convFace._face );
6938
6939   bool edgesToUpdateFound = false;
6940
6941   map< TGeomID, _EdgesOnShape* >::iterator id2eos = convFace._subIdToEOS.begin();
6942   for ( ; id2eos != convFace._subIdToEOS.end(); ++id2eos )
6943   {
6944     _EdgesOnShape& eos = * id2eos->second;
6945     if ( !eos._sWOL.IsNull() ) continue;
6946     if ( !eos._hyp.ToSmooth() ) continue;
6947     for ( size_t i = 0; i < eos._edges.size(); ++i )
6948     {
6949       _LayerEdge* ledge = eos._edges[ i ];
6950       if ( ledge->Is( _LayerEdge::UPD_NORMAL_CONV )) continue; // already checked
6951       if ( ledge->Is( _LayerEdge::MULTI_NORMAL )) continue; // not inflatable
6952
6953       gp_XYZ tgtPos = ( SMESH_NodeXYZ( ledge->_nodes[0] ) +
6954                         ledge->_normal * ledge->_lenFactor * ledge->_maxLen );
6955
6956       // the normal must be updated if distance from tgtPos to surface is less than
6957       // target thickness
6958
6959       // find an initial UV for search of a projection of tgtPos to surface
6960       const SMDS_MeshNode* nodeInFace = 0;
6961       SMDS_ElemIteratorPtr fIt = ledge->_nodes[0]->GetInverseElementIterator(SMDSAbs_Face);
6962       while ( fIt->more() && !nodeInFace )
6963       {
6964         const SMDS_MeshElement* f = fIt->next();
6965         if ( convFaceID != f->getshapeId() ) continue;
6966
6967         SMDS_ElemIteratorPtr nIt = f->nodesIterator();
6968         while ( nIt->more() && !nodeInFace )
6969         {
6970           const SMDS_MeshElement* n = nIt->next();
6971           if ( n->getshapeId() == convFaceID )
6972             nodeInFace = static_cast< const SMDS_MeshNode* >( n );
6973         }
6974       }
6975       if ( !nodeInFace )
6976         continue;
6977       gp_XY uv = helper.GetNodeUV( convFace._face, nodeInFace );
6978
6979       // projection
6980       surface->NextValueOfUV( uv, tgtPos, preci );
6981       double  dist = surface->Gap();
6982       if ( dist < 0.95 * ledge->_maxLen )
6983       {
6984         ledge->Set( _LayerEdge::UPD_NORMAL_CONV );
6985         if ( !ledge->_curvature ) ledge->_curvature = new _Curvature;
6986         ledge->_curvature->_uv.SetCoord( uv.X(), uv.Y() );
6987         edgesToUpdateFound = true;
6988       }
6989     }
6990   }
6991
6992   if ( !convFace._isTooCurved && edgesToUpdateFound )
6993   {
6994     data._convexFaces.insert( make_pair( convFaceID, convFace )).first->second;
6995   }
6996 }
6997
6998 //================================================================================
6999 /*!
7000  * \brief Modify normals of _LayerEdge's on EDGE's to avoid intersection with
7001  * _LayerEdge's on neighbor EDGE's
7002  */
7003 //================================================================================
7004
7005 bool _ViscousBuilder::updateNormals( _SolidData&         data,
7006                                      SMESH_MesherHelper& helper,
7007                                      int                 stepNb,
7008                                      double              stepSize)
7009 {
7010   updateNormalsOfC1Vertices( data );
7011
7012   if ( stepNb > 0 && !updateNormalsOfConvexFaces( data, helper, stepNb ))
7013     return false;
7014
7015   // map to store new _normal and _cosin for each intersected edge
7016   map< _LayerEdge*, _LayerEdge, _LayerEdgeCmp >           edge2newEdge;
7017   map< _LayerEdge*, _LayerEdge, _LayerEdgeCmp >::iterator e2neIt;
7018   _LayerEdge zeroEdge;
7019   zeroEdge._normal.SetCoord( 0,0,0 );
7020   zeroEdge._maxLen = Precision::Infinite();
7021   zeroEdge._nodes.resize(1); // to init _TmpMeshFaceOnEdge
7022
7023   set< _EdgesOnShape* > shapesToSmooth, edgesNoAnaSmooth;
7024
7025   double segLen, dist1, dist2, dist;
7026   vector< pair< _LayerEdge*, double > > intEdgesDist;
7027   _TmpMeshFaceOnEdge quad( &zeroEdge, &zeroEdge, 0 );
7028
7029   for ( int iter = 0; iter < 5; ++iter )
7030   {
7031     edge2newEdge.clear();
7032
7033     for ( size_t iE = 0; iE < data._collisionEdges.size(); ++iE )
7034     {
7035       _CollisionEdges& ce = data._collisionEdges[iE];
7036       _LayerEdge*   edge1 = ce._edge;
7037       if ( !edge1 /*|| edge1->Is( _LayerEdge::BLOCKED )*/) continue;
7038       _EdgesOnShape* eos1 = data.GetShapeEdges( edge1 );
7039       if ( !eos1 ) continue;
7040
7041       // detect intersections
7042       gp_Ax1 lastSeg = edge1->LastSegment( segLen, *eos1 );
7043       double testLen = 1.5 * edge1->_maxLen * edge1->_lenFactor;
7044       double     eps = 0.5;
7045       intEdgesDist.clear();
7046       double minIntDist = Precision::Infinite();
7047       for ( size_t i = 0; i < ce._intEdges.size(); i += 2 )
7048       {
7049         if ( edge1->Is( _LayerEdge::BLOCKED ) &&
7050              ce._intEdges[i  ]->Is( _LayerEdge::BLOCKED ) &&
7051              ce._intEdges[i+1]->Is( _LayerEdge::BLOCKED ))
7052           continue;
7053         double dot  = edge1->_normal * quad.GetDir( ce._intEdges[i], ce._intEdges[i+1] );
7054         double fact = ( 1.1 + dot * dot );
7055         SMESH_TNodeXYZ pSrc0( ce.nSrc(i) ), pSrc1( ce.nSrc(i+1) );
7056         SMESH_TNodeXYZ pTgt0( ce.nTgt(i) ), pTgt1( ce.nTgt(i+1) );
7057         gp_XYZ pLast0 = pSrc0 + ( pTgt0 - pSrc0 ) * fact;
7058         gp_XYZ pLast1 = pSrc1 + ( pTgt1 - pSrc1 ) * fact;
7059         dist1 = dist2 = Precision::Infinite();
7060         if ( !edge1->SegTriaInter( lastSeg, pSrc0, pLast0, pSrc1,  dist1, eps ) &&
7061              !edge1->SegTriaInter( lastSeg, pSrc1, pLast1, pLast0, dist2, eps ))
7062           continue;
7063         dist = dist1;
7064         if ( dist > testLen || dist <= 0 )
7065         {
7066           dist = dist2;
7067           if ( dist > testLen || dist <= 0 )
7068             continue;
7069         }
7070         // choose a closest edge
7071         gp_Pnt intP( lastSeg.Location().XYZ() + lastSeg.Direction().XYZ() * ( dist + segLen ));
7072         double d1 = intP.SquareDistance( pSrc0 );
7073         double d2 = intP.SquareDistance( pSrc1 );
7074         int iClose = i + ( d2 < d1 );
7075         _LayerEdge* edge2 = ce._intEdges[iClose];
7076         edge2->Unset( _LayerEdge::MARKED );
7077
7078         // choose a closest edge among neighbors
7079         gp_Pnt srcP( SMESH_TNodeXYZ( edge1->_nodes[0] ));
7080         d1 = srcP.SquareDistance( SMESH_TNodeXYZ( edge2->_nodes[0] ));
7081         for ( size_t j = 0; j < intEdgesDist.size(); ++j )
7082         {
7083           _LayerEdge * edgeJ = intEdgesDist[j].first;
7084           if ( edge2->IsNeiborOnEdge( edgeJ ))
7085           {
7086             d2 = srcP.SquareDistance( SMESH_TNodeXYZ( edgeJ->_nodes[0] ));
7087             ( d1 < d2 ? edgeJ : edge2 )->Set( _LayerEdge::MARKED );
7088           }
7089         }
7090         intEdgesDist.push_back( make_pair( edge2, dist ));
7091         // if ( Abs( d2 - d1 ) / Max( d2, d1 ) < 0.5 )
7092         // {
7093         //   iClose = i + !( d2 < d1 );
7094         //   intEdges.push_back( ce._intEdges[iClose] );
7095         //   ce._intEdges[iClose]->Unset( _LayerEdge::MARKED );
7096         // }
7097         minIntDist = Min( edge1->_len * edge1->_lenFactor - segLen + dist, minIntDist );
7098       }
7099
7100       //ce._edge = 0;
7101
7102       // compute new _normals
7103       for ( size_t i = 0; i < intEdgesDist.size(); ++i )
7104       {
7105         _LayerEdge* edge2   = intEdgesDist[i].first;
7106         double      distWgt = edge1->_len / intEdgesDist[i].second;
7107         // if ( edge1->Is( _LayerEdge::BLOCKED ) &&
7108         //      edge2->Is( _LayerEdge::BLOCKED )) continue;        
7109         if ( edge2->Is( _LayerEdge::MARKED )) continue;
7110         edge2->Set( _LayerEdge::MARKED );
7111
7112         // get a new normal
7113         gp_XYZ dir1 = edge1->_normal, dir2 = edge2->_normal;
7114
7115         double cos1 = Abs( edge1->_cosin ), cos2 = Abs( edge2->_cosin );
7116         double wgt1 = ( cos1 + 0.001 ) / ( cos1 + cos2 + 0.002 );
7117         double wgt2 = ( cos2 + 0.001 ) / ( cos1 + cos2 + 0.002 );
7118         // double cos1 = Abs( edge1->_cosin ),        cos2 = Abs( edge2->_cosin );
7119         // double sgn1 = 0.1 * ( 1 + edge1->_cosin ), sgn2 = 0.1 * ( 1 + edge2->_cosin );
7120         // double wgt1 = ( cos1 + sgn1 ) / ( cos1 + cos2 + sgn1 + sgn2 );
7121         // double wgt2 = ( cos2 + sgn2 ) / ( cos1 + cos2 + sgn1 + sgn2 );
7122         gp_XYZ newNormal = wgt1 * dir1 + wgt2 * dir2;
7123         newNormal.Normalize();
7124
7125         // get new cosin
7126         double newCos;
7127         double sgn1 = edge1->_cosin / cos1, sgn2 = edge2->_cosin / cos2;
7128         if ( cos1 < theMinSmoothCosin )
7129         {
7130           newCos = cos2 * sgn1;
7131         }
7132         else if ( cos2 > theMinSmoothCosin ) // both cos1 and cos2 > theMinSmoothCosin
7133         {
7134           newCos = ( wgt1 * cos1 + wgt2 * cos2 ) * edge1->_cosin / cos1;
7135         }
7136         else
7137         {
7138           newCos = edge1->_cosin;
7139         }
7140
7141         e2neIt = edge2newEdge.insert( make_pair( edge1, zeroEdge )).first;
7142         e2neIt->second._normal += distWgt * newNormal;
7143         e2neIt->second._cosin   = newCos;
7144         e2neIt->second.SetMaxLen( 0.7 * minIntDist / edge1->_lenFactor );
7145         if ( iter > 0 && sgn1 * sgn2 < 0 && edge1->_cosin < 0 )
7146           e2neIt->second._normal += dir2;
7147
7148         e2neIt = edge2newEdge.insert( make_pair( edge2, zeroEdge )).first;
7149         e2neIt->second._normal += distWgt * newNormal;
7150         if ( Precision::IsInfinite( zeroEdge._maxLen ))
7151         {
7152           e2neIt->second._cosin  = edge2->_cosin;
7153           e2neIt->second.SetMaxLen( 1.3 * minIntDist / edge1->_lenFactor );
7154         }
7155         if ( iter > 0 && sgn1 * sgn2 < 0 && edge2->_cosin < 0 )
7156           e2neIt->second._normal += dir1;
7157       }
7158     }
7159
7160     if ( edge2newEdge.empty() )
7161       break; //return true;
7162
7163     dumpFunction(SMESH_Comment("updateNormals")<< data._index << "_" << stepNb << "_it" << iter);
7164
7165     // Update data of edges depending on a new _normal
7166
7167     data.UnmarkEdges();
7168     for ( e2neIt = edge2newEdge.begin(); e2neIt != edge2newEdge.end(); ++e2neIt )
7169     {
7170       _LayerEdge*    edge = e2neIt->first;
7171       _LayerEdge& newEdge = e2neIt->second;
7172       _EdgesOnShape*  eos = data.GetShapeEdges( edge );
7173       if ( edge->Is( _LayerEdge::BLOCKED && newEdge._maxLen > edge->_len ))
7174         continue;
7175
7176       // Check if a new _normal is OK:
7177       newEdge._normal.Normalize();
7178       if ( !isNewNormalOk( data, *edge, newEdge._normal ))
7179       {
7180         if ( newEdge._maxLen < edge->_len && iter > 0 ) // limit _maxLen
7181         {
7182           edge->InvalidateStep( stepNb + 1, *eos, /*restoreLength=*/true  );
7183           edge->SetMaxLen( newEdge._maxLen );
7184           edge->SetNewLength( newEdge._maxLen, *eos, helper );
7185         }
7186         continue; // the new _normal is bad
7187       }
7188       // the new _normal is OK
7189
7190       // find shapes that need smoothing due to change of _normal
7191       if ( edge->_cosin   < theMinSmoothCosin &&
7192            newEdge._cosin > theMinSmoothCosin )
7193       {
7194         if ( eos->_sWOL.IsNull() )
7195         {
7196           SMDS_ElemIteratorPtr fIt = edge->_nodes[0]->GetInverseElementIterator(SMDSAbs_Face);
7197           while ( fIt->more() )
7198             shapesToSmooth.insert( data.GetShapeEdges( fIt->next()->getshapeId() ));
7199         }
7200         else // edge inflates along a FACE
7201         {
7202           TopoDS_Shape V = helper.GetSubShapeByNode( edge->_nodes[0], getMeshDS() );
7203           PShapeIteratorPtr eIt = helper.GetAncestors( V, *_mesh, TopAbs_EDGE, &eos->_sWOL );
7204           while ( const TopoDS_Shape* E = eIt->next() )
7205           {
7206             gp_Vec edgeDir = getEdgeDir( TopoDS::Edge( *E ), TopoDS::Vertex( V ));
7207             double   angle = edgeDir.Angle( newEdge._normal ); // [0,PI]
7208             if ( angle < M_PI / 2 )
7209               shapesToSmooth.insert( data.GetShapeEdges( *E ));
7210           }
7211         }
7212       }
7213
7214       double len = edge->_len;
7215       edge->InvalidateStep( stepNb + 1, *eos, /*restoreLength=*/true  );
7216       edge->SetNormal( newEdge._normal );
7217       edge->SetCosin( newEdge._cosin );
7218       edge->SetNewLength( len, *eos, helper );
7219       edge->Set( _LayerEdge::MARKED );
7220       edge->Set( _LayerEdge::NORMAL_UPDATED );
7221       edgesNoAnaSmooth.insert( eos );
7222     }
7223
7224     // Update normals and other dependent data of not intersecting _LayerEdge's
7225     // neighboring the intersecting ones
7226
7227     for ( e2neIt = edge2newEdge.begin(); e2neIt != edge2newEdge.end(); ++e2neIt )
7228     {
7229       _LayerEdge*   edge1 = e2neIt->first;
7230       _EdgesOnShape* eos1 = data.GetShapeEdges( edge1 );
7231       if ( !edge1->Is( _LayerEdge::MARKED ))
7232         continue;
7233
7234       if ( edge1->IsOnEdge() )
7235       {
7236         const SMDS_MeshNode * n1 = edge1->_2neibors->srcNode(0);
7237         const SMDS_MeshNode * n2 = edge1->_2neibors->srcNode(1);
7238         edge1->SetDataByNeighbors( n1, n2, *eos1, helper );
7239       }
7240
7241       if ( !edge1->_2neibors || !eos1->_sWOL.IsNull() )
7242         continue;
7243       for ( int j = 0; j < 2; ++j ) // loop on 2 neighbors
7244       {
7245         _LayerEdge* neighbor = edge1->_2neibors->_edges[j];
7246         if ( neighbor->Is( _LayerEdge::MARKED ) /*edge2newEdge.count ( neighbor )*/)
7247           continue; // j-th neighbor is also intersected
7248         _LayerEdge* prevEdge = edge1;
7249         const int nbSteps = 10;
7250         for ( int step = nbSteps; step; --step ) // step from edge1 in j-th direction
7251         {
7252           if ( neighbor->Is( _LayerEdge::BLOCKED ) ||
7253                neighbor->Is( _LayerEdge::MARKED ))
7254             break;
7255           _EdgesOnShape* eos = data.GetShapeEdges( neighbor );
7256           if ( !eos ) continue;
7257           _LayerEdge* nextEdge = neighbor;
7258           if ( neighbor->_2neibors )
7259           {
7260             int iNext = 0;
7261             nextEdge = neighbor->_2neibors->_edges[iNext];
7262             if ( nextEdge == prevEdge )
7263               nextEdge = neighbor->_2neibors->_edges[ ++iNext ];
7264           }
7265           double r = double(step-1)/nbSteps/(iter+1);
7266           if ( !nextEdge->_2neibors )
7267             r = Min( r, 0.5 );
7268
7269           gp_XYZ newNorm = prevEdge->_normal * r + nextEdge->_normal * (1-r);
7270           newNorm.Normalize();
7271           if ( !isNewNormalOk( data, *neighbor, newNorm ))
7272             break;
7273
7274           double len = neighbor->_len;
7275           neighbor->InvalidateStep( stepNb + 1, *eos, /*restoreLength=*/true  );
7276           neighbor->SetNormal( newNorm );
7277           neighbor->SetCosin( prevEdge->_cosin * r + nextEdge->_cosin * (1-r) );
7278           if ( neighbor->_2neibors )
7279             neighbor->SetDataByNeighbors( prevEdge->_nodes[0], nextEdge->_nodes[0], *eos, helper );
7280           neighbor->SetNewLength( len, *eos, helper );
7281           neighbor->Set( _LayerEdge::MARKED );
7282           neighbor->Set( _LayerEdge::NORMAL_UPDATED );
7283           edgesNoAnaSmooth.insert( eos );
7284
7285           if ( !neighbor->_2neibors )
7286             break; // neighbor is on VERTEX
7287
7288           // goto the next neighbor
7289           prevEdge = neighbor;
7290           neighbor = nextEdge;
7291         }
7292       }
7293     }
7294     dumpFunctionEnd();
7295   } // iterations
7296
7297   data.AddShapesToSmooth( shapesToSmooth, &edgesNoAnaSmooth );
7298
7299   return true;
7300 }
7301
7302 //================================================================================
7303 /*!
7304  * \brief Check if a new normal is OK
7305  */
7306 //================================================================================
7307
7308 bool _ViscousBuilder::isNewNormalOk( _SolidData&   data,
7309                                      _LayerEdge&   edge,
7310                                      const gp_XYZ& newNormal)
7311 {
7312   // check a min angle between the newNormal and surrounding faces
7313   vector<_Simplex> simplices;
7314   SMESH_TNodeXYZ n0( edge._nodes[0] ), n1, n2;
7315   _Simplex::GetSimplices( n0._node, simplices, data._ignoreFaceIds, &data );
7316   double newMinDot = 1, curMinDot = 1;
7317   for ( size_t i = 0; i < simplices.size(); ++i )
7318   {
7319     n1.Set( simplices[i]._nPrev );
7320     n2.Set( simplices[i]._nNext );
7321     gp_XYZ normFace = ( n1 - n0 ) ^ ( n2 - n0 );
7322     double normLen2 = normFace.SquareModulus();
7323     if ( normLen2 < std::numeric_limits<double>::min() )
7324       continue;
7325     normFace /= Sqrt( normLen2 );
7326     newMinDot = Min( newNormal    * normFace, newMinDot );
7327     curMinDot = Min( edge._normal * normFace, curMinDot );
7328   }
7329   bool ok = true;
7330   if ( newMinDot < 0.5 )
7331   {
7332     ok = ( newMinDot >= curMinDot * 0.9 );
7333     //return ( newMinDot >= ( curMinDot * ( 0.8 + 0.1 * edge.NbSteps() )));
7334     // double initMinDot2 = 1. - edge._cosin * edge._cosin;
7335     // return ( newMinDot * newMinDot ) >= ( 0.8 * initMinDot2 );
7336   }
7337
7338   return ok;
7339 }
7340
7341 //================================================================================
7342 /*!
7343  * \brief Modify normals of _LayerEdge's on FACE to reflex smoothing
7344  */
7345 //================================================================================
7346
7347 bool _ViscousBuilder::updateNormalsOfSmoothed( _SolidData&         data,
7348                                                SMESH_MesherHelper& helper,
7349                                                const int           nbSteps,
7350                                                const double        stepSize )
7351 {
7352   if ( data._nbShapesToSmooth == 0 || nbSteps == 0 )
7353     return true; // no shapes needing smoothing
7354
7355   for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
7356   {
7357     _EdgesOnShape& eos = data._edgesOnShape[ iS ];
7358     if ( //!eos._toSmooth ||  _eosC1 have _toSmooth == false
7359          !eos._hyp.ToSmooth() ||
7360          eos.ShapeType() != TopAbs_FACE ||
7361          eos._edges.empty() )
7362       continue;
7363
7364     bool toSmooth = ( eos._edges[ 0 ]->NbSteps() >= nbSteps+1 );
7365     if ( !toSmooth ) continue;
7366
7367     for ( size_t i = 0; i < eos._edges.size(); ++i )
7368     {
7369       _LayerEdge* edge = eos._edges[i];
7370       if ( !edge->Is( _LayerEdge::SMOOTHED ))
7371         continue;
7372       if ( edge->Is( _LayerEdge::DIFFICULT ) && nbSteps != 1 )
7373         continue;
7374
7375       const gp_XYZ& pPrev = edge->PrevPos();
7376       const gp_XYZ& pLast = edge->_pos.back();
7377       gp_XYZ      stepVec = pLast - pPrev;
7378       double realStepSize = stepVec.Modulus();
7379       if ( realStepSize < numeric_limits<double>::min() )
7380         continue;
7381
7382       edge->_lenFactor = realStepSize / stepSize;
7383       edge->_normal    = stepVec / realStepSize;
7384       edge->Set( _LayerEdge::NORMAL_UPDATED );
7385     }
7386   }
7387
7388   return true;
7389 }
7390
7391 //================================================================================
7392 /*!
7393  * \brief Modify normals of _LayerEdge's on C1 VERTEXes
7394  */
7395 //================================================================================
7396
7397 void _ViscousBuilder::updateNormalsOfC1Vertices( _SolidData& data )
7398 {
7399   for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
7400   {
7401     _EdgesOnShape& eov = data._edgesOnShape[ iS ];
7402     if ( eov._eosC1.empty() ||
7403          eov.ShapeType() != TopAbs_VERTEX ||
7404          eov._edges.empty() )
7405       continue;
7406
7407     gp_XYZ newNorm   = eov._edges[0]->_normal;
7408     double curThick  = eov._edges[0]->_len * eov._edges[0]->_lenFactor;
7409     bool normChanged = false;
7410
7411     for ( size_t i = 0; i < eov._eosC1.size(); ++i )
7412     {
7413       _EdgesOnShape*   eoe = eov._eosC1[i];
7414       const TopoDS_Edge& e = TopoDS::Edge( eoe->_shape );
7415       const double    eLen = SMESH_Algo::EdgeLength( e );
7416       TopoDS_Shape    oppV = SMESH_MesherHelper::IthVertex( 0, e );
7417       if ( oppV.IsSame( eov._shape ))
7418         oppV = SMESH_MesherHelper::IthVertex( 1, e );
7419       _EdgesOnShape* eovOpp = data.GetShapeEdges( oppV );
7420       if ( !eovOpp || eovOpp->_edges.empty() ) continue;
7421       if ( eov._edges[0]->Is( _LayerEdge::BLOCKED )) continue;
7422
7423       double curThickOpp = eovOpp->_edges[0]->_len * eovOpp->_edges[0]->_lenFactor;
7424       if ( curThickOpp + curThick < eLen )
7425         continue;
7426
7427       double wgt = 2. * curThick / eLen;
7428       newNorm += wgt * eovOpp->_edges[0]->_normal;
7429       normChanged = true;
7430     }
7431     if ( normChanged )
7432     {
7433       eov._edges[0]->SetNormal( newNorm.Normalized() );
7434       eov._edges[0]->Set( _LayerEdge::NORMAL_UPDATED );
7435     }
7436   }
7437 }
7438
7439 //================================================================================
7440 /*!
7441  * \brief Modify normals of _LayerEdge's on _ConvexFace's
7442  */
7443 //================================================================================
7444
7445 bool _ViscousBuilder::updateNormalsOfConvexFaces( _SolidData&         data,
7446                                                   SMESH_MesherHelper& helper,
7447                                                   int                 stepNb )
7448 {
7449   SMESHDS_Mesh* meshDS = helper.GetMeshDS();
7450   bool isOK;
7451
7452   map< TGeomID, _ConvexFace >::iterator id2face = data._convexFaces.begin();
7453   for ( ; id2face != data._convexFaces.end(); ++id2face )
7454   {
7455     _ConvexFace & convFace = (*id2face).second;
7456     convFace._normalsFixedOnBorders = false; // to update at each inflation step
7457
7458     if ( convFace._normalsFixed )
7459       continue; // already fixed
7460     if ( convFace.CheckPrisms() )
7461       continue; // nothing to fix
7462
7463     convFace._normalsFixed = true;
7464
7465     BRepAdaptor_Surface surface ( convFace._face, false );
7466     BRepLProp_SLProps   surfProp( surface, 2, 1e-6 );
7467
7468     // check if the convex FACE is of spherical shape
7469
7470     Bnd_B3d centersBox; // bbox of centers of curvature of _LayerEdge's on VERTEXes
7471     Bnd_B3d nodesBox;
7472     gp_Pnt  center;
7473
7474     map< TGeomID, _EdgesOnShape* >::iterator id2eos = convFace._subIdToEOS.begin();
7475     for ( ; id2eos != convFace._subIdToEOS.end(); ++id2eos )
7476     {
7477       _EdgesOnShape& eos = *(id2eos->second);
7478       if ( eos.ShapeType() == TopAbs_VERTEX )
7479       {
7480         _LayerEdge* ledge = eos._edges[ 0 ];
7481         if ( convFace.GetCenterOfCurvature( ledge, surfProp, helper, center ))
7482           centersBox.Add( center );
7483       }
7484       for ( size_t i = 0; i < eos._edges.size(); ++i )
7485         nodesBox.Add( SMESH_TNodeXYZ( eos._edges[ i ]->_nodes[0] ));
7486     }
7487     if ( centersBox.IsVoid() )
7488     {
7489       debugMsg( "Error: centersBox.IsVoid()" );
7490       return false;
7491     }
7492     const bool isSpherical =
7493       ( centersBox.SquareExtent() < 1e-6 * nodesBox.SquareExtent() );
7494
7495     int nbEdges = helper.Count( convFace._face, TopAbs_EDGE, /*ignoreSame=*/false );
7496     vector < _CentralCurveOnEdge > centerCurves( nbEdges );
7497
7498     if ( isSpherical )
7499     {
7500       // set _LayerEdge::_normal as average of all normals
7501
7502       // WARNING: different density of nodes on EDGEs is not taken into account that
7503       // can lead to an improper new normal
7504
7505       gp_XYZ avgNormal( 0,0,0 );
7506       nbEdges = 0;
7507       id2eos = convFace._subIdToEOS.begin();
7508       for ( ; id2eos != convFace._subIdToEOS.end(); ++id2eos )
7509       {
7510         _EdgesOnShape& eos = *(id2eos->second);
7511         // set data of _CentralCurveOnEdge
7512         if ( eos.ShapeType() == TopAbs_EDGE )
7513         {
7514           _CentralCurveOnEdge& ceCurve = centerCurves[ nbEdges++ ];
7515           ceCurve.SetShapes( TopoDS::Edge( eos._shape ), convFace, data, helper );
7516           if ( !eos._sWOL.IsNull() )
7517             ceCurve._adjFace.Nullify();
7518           else
7519             ceCurve._ledges.insert( ceCurve._ledges.end(),
7520                                     eos._edges.begin(), eos._edges.end());
7521         }
7522         // summarize normals
7523         for ( size_t i = 0; i < eos._edges.size(); ++i )
7524           avgNormal += eos._edges[ i ]->_normal;
7525       }
7526       double normSize = avgNormal.SquareModulus();
7527       if ( normSize < 1e-200 )
7528       {
7529         debugMsg( "updateNormalsOfConvexFaces(): zero avgNormal" );
7530         return false;
7531       }
7532       avgNormal /= Sqrt( normSize );
7533
7534       // compute new _LayerEdge::_cosin on EDGEs
7535       double avgCosin = 0;
7536       int     nbCosin = 0;
7537       gp_Vec inFaceDir;
7538       for ( size_t iE = 0; iE < centerCurves.size(); ++iE )
7539       {
7540         _CentralCurveOnEdge& ceCurve = centerCurves[ iE ];
7541         if ( ceCurve._adjFace.IsNull() )
7542           continue;
7543         for ( size_t iLE = 0; iLE < ceCurve._ledges.size(); ++iLE )
7544         {
7545           const SMDS_MeshNode* node = ceCurve._ledges[ iLE ]->_nodes[0];
7546           inFaceDir = getFaceDir( ceCurve._adjFace, ceCurve._edge, node, helper, isOK );
7547           if ( isOK )
7548           {
7549             double angle = inFaceDir.Angle( avgNormal ); // [0,PI]
7550             ceCurve._ledges[ iLE ]->_cosin = Cos( angle );
7551             avgCosin += ceCurve._ledges[ iLE ]->_cosin;
7552             nbCosin++;
7553           }
7554         }
7555       }
7556       if ( nbCosin > 0 )
7557         avgCosin /= nbCosin;
7558
7559       // set _LayerEdge::_normal = avgNormal
7560       id2eos = convFace._subIdToEOS.begin();
7561       for ( ; id2eos != convFace._subIdToEOS.end(); ++id2eos )
7562       {
7563         _EdgesOnShape& eos = *(id2eos->second);
7564         if ( eos.ShapeType() != TopAbs_EDGE )
7565           for ( size_t i = 0; i < eos._edges.size(); ++i )
7566             eos._edges[ i ]->_cosin = avgCosin;
7567
7568         for ( size_t i = 0; i < eos._edges.size(); ++i )
7569         {
7570           eos._edges[ i ]->SetNormal( avgNormal );
7571           eos._edges[ i ]->Set( _LayerEdge::NORMAL_UPDATED );
7572         }
7573       }
7574     }
7575     else // if ( isSpherical )
7576     {
7577       // We suppose that centers of curvature at all points of the FACE
7578       // lie on some curve, let's call it "central curve". For all _LayerEdge's
7579       // having a common center of curvature we define the same new normal
7580       // as a sum of normals of _LayerEdge's on EDGEs among them.
7581
7582       // get all centers of curvature for each EDGE
7583
7584       helper.SetSubShape( convFace._face );
7585       _LayerEdge* vertexLEdges[2], **edgeLEdge, **edgeLEdgeEnd;
7586
7587       TopExp_Explorer edgeExp( convFace._face, TopAbs_EDGE );
7588       for ( int iE = 0; edgeExp.More(); edgeExp.Next(), ++iE )
7589       {
7590         const TopoDS_Edge& edge = TopoDS::Edge( edgeExp.Current() );
7591
7592         // set adjacent FACE
7593         centerCurves[ iE ].SetShapes( edge, convFace, data, helper );
7594
7595         // get _LayerEdge's of the EDGE
7596         TGeomID edgeID = meshDS->ShapeToIndex( edge );
7597         _EdgesOnShape* eos = data.GetShapeEdges( edgeID );
7598         if ( !eos || eos->_edges.empty() )
7599         {
7600           // no _LayerEdge's on EDGE, use _LayerEdge's on VERTEXes
7601           for ( int iV = 0; iV < 2; ++iV )
7602           {
7603             TopoDS_Vertex v = helper.IthVertex( iV, edge );
7604             TGeomID     vID = meshDS->ShapeToIndex( v );
7605             eos = data.GetShapeEdges( vID );
7606             vertexLEdges[ iV ] = eos->_edges[ 0 ];
7607           }
7608           edgeLEdge    = &vertexLEdges[0];
7609           edgeLEdgeEnd = edgeLEdge + 2;
7610
7611           centerCurves[ iE ]._adjFace.Nullify();
7612         }
7613         else
7614         {
7615           if ( ! eos->_toSmooth )
7616             data.SortOnEdge( edge, eos->_edges );
7617           edgeLEdge    = &eos->_edges[ 0 ];
7618           edgeLEdgeEnd = edgeLEdge + eos->_edges.size();
7619           vertexLEdges[0] = eos->_edges.front()->_2neibors->_edges[0];
7620           vertexLEdges[1] = eos->_edges.back() ->_2neibors->_edges[1];
7621
7622           if ( ! eos->_sWOL.IsNull() )
7623             centerCurves[ iE ]._adjFace.Nullify();
7624         }
7625
7626         // Get curvature centers
7627
7628         centersBox.Clear();
7629
7630         if ( edgeLEdge[0]->IsOnEdge() &&
7631              convFace.GetCenterOfCurvature( vertexLEdges[0], surfProp, helper, center ))
7632         { // 1st VERTEX
7633           centerCurves[ iE ].Append( center, vertexLEdges[0] );
7634           centersBox.Add( center );
7635         }
7636         for ( ; edgeLEdge < edgeLEdgeEnd; ++edgeLEdge )
7637           if ( convFace.GetCenterOfCurvature( *edgeLEdge, surfProp, helper, center ))
7638           { // EDGE or VERTEXes
7639             centerCurves[ iE ].Append( center, *edgeLEdge );
7640             centersBox.Add( center );
7641           }
7642         if ( edgeLEdge[-1]->IsOnEdge() &&
7643              convFace.GetCenterOfCurvature( vertexLEdges[1], surfProp, helper, center ))
7644         { // 2nd VERTEX
7645           centerCurves[ iE ].Append( center, vertexLEdges[1] );
7646           centersBox.Add( center );
7647         }
7648         centerCurves[ iE ]._isDegenerated =
7649           ( centersBox.IsVoid() || centersBox.SquareExtent() < 1e-6 * nodesBox.SquareExtent() );
7650
7651       } // loop on EDGES of convFace._face to set up data of centerCurves
7652
7653       // Compute new normals for _LayerEdge's on EDGEs
7654
7655       double avgCosin = 0;
7656       int     nbCosin = 0;
7657       gp_Vec inFaceDir;
7658       for ( size_t iE1 = 0; iE1 < centerCurves.size(); ++iE1 )
7659       {
7660         _CentralCurveOnEdge& ceCurve = centerCurves[ iE1 ];
7661         if ( ceCurve._isDegenerated )
7662           continue;
7663         const vector< gp_Pnt >& centers = ceCurve._curvaCenters;
7664         vector< gp_XYZ > &   newNormals = ceCurve._normals;
7665         for ( size_t iC1 = 0; iC1 < centers.size(); ++iC1 )
7666         {
7667           isOK = false;
7668           for ( size_t iE2 = 0; iE2 < centerCurves.size() && !isOK; ++iE2 )
7669           {
7670             if ( iE1 != iE2 )
7671               isOK = centerCurves[ iE2 ].FindNewNormal( centers[ iC1 ], newNormals[ iC1 ]);
7672           }
7673           if ( isOK && !ceCurve._adjFace.IsNull() )
7674           {
7675             // compute new _LayerEdge::_cosin
7676             const SMDS_MeshNode* node = ceCurve._ledges[ iC1 ]->_nodes[0];
7677             inFaceDir = getFaceDir( ceCurve._adjFace, ceCurve._edge, node, helper, isOK );
7678             if ( isOK )
7679             {
7680               double angle = inFaceDir.Angle( newNormals[ iC1 ] ); // [0,PI]
7681               ceCurve._ledges[ iC1 ]->_cosin = Cos( angle );
7682               avgCosin += ceCurve._ledges[ iC1 ]->_cosin;
7683               nbCosin++;
7684             }
7685           }
7686         }
7687       }
7688       // set new normals to _LayerEdge's of NOT degenerated central curves
7689       for ( size_t iE = 0; iE < centerCurves.size(); ++iE )
7690       {
7691         if ( centerCurves[ iE ]._isDegenerated )
7692           continue;
7693         for ( size_t iLE = 0; iLE < centerCurves[ iE ]._ledges.size(); ++iLE )
7694         {
7695           centerCurves[ iE ]._ledges[ iLE ]->SetNormal( centerCurves[ iE ]._normals[ iLE ]);
7696           centerCurves[ iE ]._ledges[ iLE ]->Set( _LayerEdge::NORMAL_UPDATED );
7697         }
7698       }
7699       // set new normals to _LayerEdge's of     degenerated central curves
7700       for ( size_t iE = 0; iE < centerCurves.size(); ++iE )
7701       {
7702         if ( !centerCurves[ iE ]._isDegenerated ||
7703              centerCurves[ iE ]._ledges.size() < 3 )
7704           continue;
7705         // new normal is an average of new normals at VERTEXes that
7706         // was computed on non-degenerated _CentralCurveOnEdge's
7707         gp_XYZ newNorm = ( centerCurves[ iE ]._ledges.front()->_normal +
7708                            centerCurves[ iE ]._ledges.back ()->_normal );
7709         double sz = newNorm.Modulus();
7710         if ( sz < 1e-200 )
7711           continue;
7712         newNorm /= sz;
7713         double newCosin = ( 0.5 * centerCurves[ iE ]._ledges.front()->_cosin +
7714                             0.5 * centerCurves[ iE ]._ledges.back ()->_cosin );
7715         for ( size_t iLE = 1, nb = centerCurves[ iE ]._ledges.size() - 1; iLE < nb; ++iLE )
7716         {
7717           centerCurves[ iE ]._ledges[ iLE ]->SetNormal( newNorm );
7718           centerCurves[ iE ]._ledges[ iLE ]->_cosin   = newCosin;
7719           centerCurves[ iE ]._ledges[ iLE ]->Set( _LayerEdge::NORMAL_UPDATED );
7720         }
7721       }
7722
7723       // Find new normals for _LayerEdge's based on FACE
7724
7725       if ( nbCosin > 0 )
7726         avgCosin /= nbCosin;
7727       const TGeomID faceID = meshDS->ShapeToIndex( convFace._face );
7728       map< TGeomID, _EdgesOnShape* >::iterator id2eos = convFace._subIdToEOS.find( faceID );
7729       if ( id2eos != convFace._subIdToEOS.end() )
7730       {
7731         int iE = 0;
7732         gp_XYZ newNorm;
7733         _EdgesOnShape& eos = * ( id2eos->second );
7734         for ( size_t i = 0; i < eos._edges.size(); ++i )
7735         {
7736           _LayerEdge* ledge = eos._edges[ i ];
7737           if ( !convFace.GetCenterOfCurvature( ledge, surfProp, helper, center ))
7738             continue;
7739           for ( size_t i = 0; i < centerCurves.size(); ++i, ++iE )
7740           {
7741             iE = iE % centerCurves.size();
7742             if ( centerCurves[ iE ]._isDegenerated )
7743               continue;
7744             newNorm.SetCoord( 0,0,0 );
7745             if ( centerCurves[ iE ].FindNewNormal( center, newNorm ))
7746             {
7747               ledge->SetNormal( newNorm );
7748               ledge->_cosin  = avgCosin;
7749               ledge->Set( _LayerEdge::NORMAL_UPDATED );
7750               break;
7751             }
7752           }
7753         }
7754       }
7755
7756     } // not a quasi-spherical FACE
7757
7758     // Update _LayerEdge's data according to a new normal
7759
7760     dumpFunction(SMESH_Comment("updateNormalsOfConvexFaces")<<data._index
7761                  <<"_F"<<meshDS->ShapeToIndex( convFace._face ));
7762
7763     id2eos = convFace._subIdToEOS.begin();
7764     for ( ; id2eos != convFace._subIdToEOS.end(); ++id2eos )
7765     {
7766       _EdgesOnShape& eos = * ( id2eos->second );
7767       for ( size_t i = 0; i < eos._edges.size(); ++i )
7768       {
7769         _LayerEdge* & ledge = eos._edges[ i ];
7770         double len = ledge->_len;
7771         ledge->InvalidateStep( stepNb + 1, eos, /*restoreLength=*/true );
7772         ledge->SetCosin( ledge->_cosin );
7773         ledge->SetNewLength( len, eos, helper );
7774       }
7775       if ( eos.ShapeType() != TopAbs_FACE )
7776         for ( size_t i = 0; i < eos._edges.size(); ++i )
7777         {
7778           _LayerEdge* ledge = eos._edges[ i ];
7779           for ( size_t iN = 0; iN < ledge->_neibors.size(); ++iN )
7780           {
7781             _LayerEdge* neibor = ledge->_neibors[iN];
7782             if ( neibor->_nodes[0]->GetPosition()->GetDim() == 2 )
7783             {
7784               neibor->Set( _LayerEdge::NEAR_BOUNDARY );
7785               neibor->Set( _LayerEdge::MOVED );
7786               neibor->SetSmooLen( neibor->_len );
7787             }
7788           }
7789         }
7790     } // loop on sub-shapes of convFace._face
7791
7792     // Find FACEs adjacent to convFace._face that got necessity to smooth
7793     // as a result of normals modification
7794
7795     set< _EdgesOnShape* > adjFacesToSmooth;
7796     for ( size_t iE = 0; iE < centerCurves.size(); ++iE )
7797     {
7798       if ( centerCurves[ iE ]._adjFace.IsNull() ||
7799            centerCurves[ iE ]._adjFaceToSmooth )
7800         continue;
7801       for ( size_t iLE = 0; iLE < centerCurves[ iE ]._ledges.size(); ++iLE )
7802       {
7803         if ( centerCurves[ iE ]._ledges[ iLE ]->_cosin > theMinSmoothCosin )
7804         {
7805           adjFacesToSmooth.insert( data.GetShapeEdges( centerCurves[ iE ]._adjFace ));
7806           break;
7807         }
7808       }
7809     }
7810     data.AddShapesToSmooth( adjFacesToSmooth );
7811
7812     dumpFunctionEnd();
7813
7814
7815   } // loop on data._convexFaces
7816
7817   return true;
7818 }
7819
7820 //================================================================================
7821 /*!
7822  * \brief Return max curvature of a FACE
7823  */
7824 //================================================================================
7825
7826 double _ConvexFace::GetMaxCurvature( _SolidData&         data,
7827                                      _EdgesOnShape&      eof,
7828                                      BRepLProp_SLProps&  surfProp,
7829                                      SMESH_MesherHelper& helper)
7830 {
7831   double maxCurvature = 0;
7832
7833   TopoDS_Face F = TopoDS::Face( eof._shape );
7834
7835   const int           nbTestPnt = 5;
7836   const double        oriFactor = ( F.Orientation() == TopAbs_REVERSED ? +1. : -1. );
7837   SMESH_subMeshIteratorPtr smIt = eof._subMesh->getDependsOnIterator(/*includeSelf=*/true);
7838   while ( smIt->more() )
7839   {
7840     SMESH_subMesh* sm = smIt->next();
7841     const TGeomID subID = sm->GetId();
7842
7843     // find _LayerEdge's of a sub-shape
7844     _EdgesOnShape* eos;
7845     if (( eos = data.GetShapeEdges( subID )))
7846       this->_subIdToEOS.insert( make_pair( subID, eos ));
7847     else
7848       continue;
7849
7850     // check concavity and curvature and limit data._stepSize
7851     const double minCurvature =
7852       1. / ( eos->_hyp.GetTotalThickness() * ( 1 + theThickToIntersection ));
7853     size_t iStep = Max( 1, eos->_edges.size() / nbTestPnt );
7854     for ( size_t i = 0; i < eos->_edges.size(); i += iStep )
7855     {
7856       gp_XY uv = helper.GetNodeUV( F, eos->_edges[ i ]->_nodes[0] );
7857       surfProp.SetParameters( uv.X(), uv.Y() );
7858       if ( surfProp.IsCurvatureDefined() )
7859       {
7860         double curvature = Max( surfProp.MaxCurvature() * oriFactor,
7861                                 surfProp.MinCurvature() * oriFactor );
7862         maxCurvature = Max( maxCurvature, curvature );
7863
7864         if ( curvature > minCurvature )
7865           this->_isTooCurved = true;
7866       }
7867     }
7868   } // loop on sub-shapes of the FACE
7869
7870   return maxCurvature;
7871 }
7872
7873 //================================================================================
7874 /*!
7875  * \brief Finds a center of curvature of a surface at a _LayerEdge
7876  */
7877 //================================================================================
7878
7879 bool _ConvexFace::GetCenterOfCurvature( _LayerEdge*         ledge,
7880                                         BRepLProp_SLProps&  surfProp,
7881                                         SMESH_MesherHelper& helper,
7882                                         gp_Pnt &            center ) const
7883 {
7884   gp_XY uv = helper.GetNodeUV( _face, ledge->_nodes[0] );
7885   surfProp.SetParameters( uv.X(), uv.Y() );
7886   if ( !surfProp.IsCurvatureDefined() )
7887     return false;
7888
7889   const double oriFactor = ( _face.Orientation() == TopAbs_REVERSED ? +1. : -1. );
7890   double surfCurvatureMax = surfProp.MaxCurvature() * oriFactor;
7891   double surfCurvatureMin = surfProp.MinCurvature() * oriFactor;
7892   if ( surfCurvatureMin > surfCurvatureMax )
7893     center = surfProp.Value().Translated( surfProp.Normal().XYZ() / surfCurvatureMin * oriFactor );
7894   else
7895     center = surfProp.Value().Translated( surfProp.Normal().XYZ() / surfCurvatureMax * oriFactor );
7896
7897   return true;
7898 }
7899
7900 //================================================================================
7901 /*!
7902  * \brief Check that prisms are not distorted
7903  */
7904 //================================================================================
7905
7906 bool _ConvexFace::CheckPrisms() const
7907 {
7908   double vol = 0;
7909   for ( size_t i = 0; i < _simplexTestEdges.size(); ++i )
7910   {
7911     const _LayerEdge* edge = _simplexTestEdges[i];
7912     SMESH_TNodeXYZ tgtXYZ( edge->_nodes.back() );
7913     for ( size_t j = 0; j < edge->_simplices.size(); ++j )
7914       if ( !edge->_simplices[j].IsForward( edge->_nodes[0], tgtXYZ, vol ))
7915       {
7916         debugMsg( "Bad simplex of _simplexTestEdges ("
7917                   << " "<< edge->_nodes[0]->GetID()<< " "<< tgtXYZ._node->GetID()
7918                   << " "<< edge->_simplices[j]._nPrev->GetID()
7919                   << " "<< edge->_simplices[j]._nNext->GetID() << " )" );
7920         return false;
7921       }
7922   }
7923   return true;
7924 }
7925
7926 //================================================================================
7927 /*!
7928  * \brief Try to compute a new normal by interpolating normals of _LayerEdge's
7929  *        stored in this _CentralCurveOnEdge.
7930  *  \param [in] center - curvature center of a point of another _CentralCurveOnEdge.
7931  *  \param [in,out] newNormal - current normal at this point, to be redefined
7932  *  \return bool - true if succeeded.
7933  */
7934 //================================================================================
7935
7936 bool _CentralCurveOnEdge::FindNewNormal( const gp_Pnt& center, gp_XYZ& newNormal )
7937 {
7938   if ( this->_isDegenerated )
7939     return false;
7940
7941   // find two centers the given one lies between
7942
7943   for ( size_t i = 0, nb = _curvaCenters.size()-1;  i < nb;  ++i )
7944   {
7945     double sl2 = 1.001 * _segLength2[ i ];
7946
7947     double d1 = center.SquareDistance( _curvaCenters[ i ]);
7948     if ( d1 > sl2 )
7949       continue;
7950     
7951     double d2 = center.SquareDistance( _curvaCenters[ i+1 ]);
7952     if ( d2 > sl2 || d2 + d1 < 1e-100 )
7953       continue;
7954
7955     d1 = Sqrt( d1 );
7956     d2 = Sqrt( d2 );
7957     double r = d1 / ( d1 + d2 );
7958     gp_XYZ norm = (( 1. - r ) * _ledges[ i   ]->_normal +
7959                    (      r ) * _ledges[ i+1 ]->_normal );
7960     norm.Normalize();
7961
7962     newNormal += norm;
7963     double sz = newNormal.Modulus();
7964     if ( sz < 1e-200 )
7965       break;
7966     newNormal /= sz;
7967     return true;
7968   }
7969   return false;
7970 }
7971
7972 //================================================================================
7973 /*!
7974  * \brief Set shape members
7975  */
7976 //================================================================================
7977
7978 void _CentralCurveOnEdge::SetShapes( const TopoDS_Edge&  edge,
7979                                      const _ConvexFace&  convFace,
7980                                      _SolidData&         data,
7981                                      SMESH_MesherHelper& helper)
7982 {
7983   _edge = edge;
7984
7985   PShapeIteratorPtr fIt = helper.GetAncestors( edge, *helper.GetMesh(), TopAbs_FACE );
7986   while ( const TopoDS_Shape* F = fIt->next())
7987     if ( !convFace._face.IsSame( *F ))
7988     {
7989       _adjFace = TopoDS::Face( *F );
7990       _adjFaceToSmooth = false;
7991       // _adjFace already in a smoothing queue ?
7992       if ( _EdgesOnShape* eos = data.GetShapeEdges( _adjFace ))
7993         _adjFaceToSmooth = eos->_toSmooth;
7994       break;
7995     }
7996 }
7997
7998 //================================================================================
7999 /*!
8000  * \brief Looks for intersection of it's last segment with faces
8001  *  \param distance - returns shortest distance from the last node to intersection
8002  */
8003 //================================================================================
8004
8005 bool _LayerEdge::FindIntersection( SMESH_ElementSearcher&   searcher,
8006                                    double &                 distance,
8007                                    const double&            epsilon,
8008                                    _EdgesOnShape&           eos,
8009                                    const SMDS_MeshElement** intFace)
8010 {
8011   vector< const SMDS_MeshElement* > suspectFaces;
8012   double segLen;
8013   gp_Ax1 lastSegment = LastSegment( segLen, eos );
8014   searcher.GetElementsNearLine( lastSegment, SMDSAbs_Face, suspectFaces );
8015
8016   bool segmentIntersected = false;
8017   distance = Precision::Infinite();
8018   int iFace = -1; // intersected face
8019   for ( size_t j = 0 ; j < suspectFaces.size() /*&& !segmentIntersected*/; ++j )
8020   {
8021     const SMDS_MeshElement* face = suspectFaces[j];
8022     if ( face->GetNodeIndex( _nodes.back() ) >= 0 ||
8023          face->GetNodeIndex( _nodes[0]     ) >= 0 )
8024       continue; // face sharing _LayerEdge node
8025     const int nbNodes = face->NbCornerNodes();
8026     bool intFound = false;
8027     double dist;
8028     SMDS_MeshElement::iterator nIt = face->begin_nodes();
8029     if ( nbNodes == 3 )
8030     {
8031       intFound = SegTriaInter( lastSegment, *nIt++, *nIt++, *nIt++, dist, epsilon );
8032     }
8033     else
8034     {
8035       const SMDS_MeshNode* tria[3];
8036       tria[0] = *nIt++;
8037       tria[1] = *nIt++;
8038       for ( int n2 = 2; n2 < nbNodes && !intFound; ++n2 )
8039       {
8040         tria[2] = *nIt++;
8041         intFound = SegTriaInter(lastSegment, tria[0], tria[1], tria[2], dist, epsilon );
8042         tria[1] = tria[2];
8043       }
8044     }
8045     if ( intFound )
8046     {
8047       if ( dist < segLen*(1.01) && dist > -(_len*_lenFactor-segLen) )
8048         segmentIntersected = true;
8049       if ( distance > dist )
8050         distance = dist, iFace = j;
8051     }
8052   }
8053   if ( intFace ) *intFace = ( iFace != -1 ) ? suspectFaces[iFace] : 0;
8054
8055   distance -= segLen;
8056
8057   if ( segmentIntersected )
8058   {
8059 #ifdef __myDEBUG
8060     SMDS_MeshElement::iterator nIt = suspectFaces[iFace]->begin_nodes();
8061     gp_XYZ intP( lastSegment.Location().XYZ() + lastSegment.Direction().XYZ() * ( distance+segLen ));
8062     cout << "nodes: tgt " << _nodes.back()->GetID() << " src " << _nodes[0]->GetID()
8063          << ", intersection with face ("
8064          << (*nIt++)->GetID()<<" "<< (*nIt++)->GetID()<<" "<< (*nIt++)->GetID()
8065          << ") at point (" << intP.X() << ", " << intP.Y() << ", " << intP.Z()
8066          << ") distance = " << distance << endl;
8067 #endif
8068   }
8069
8070   return segmentIntersected;
8071 }
8072
8073 //================================================================================
8074 /*!
8075  * \brief Returns a point used to check orientation of _simplices
8076  */
8077 //================================================================================
8078
8079 gp_XYZ _LayerEdge::PrevCheckPos( _EdgesOnShape* eos ) const
8080 {
8081   size_t i = Is( NORMAL_UPDATED ) && IsOnFace() ? _pos.size()-2 : 0;
8082
8083   if ( !eos || eos->_sWOL.IsNull() )
8084     return _pos[ i ];
8085
8086   if ( eos->SWOLType() == TopAbs_EDGE )
8087   {
8088     return BRepAdaptor_Curve( TopoDS::Edge( eos->_sWOL )).Value( _pos[i].X() ).XYZ();
8089   }
8090   //else //  TopAbs_FACE
8091
8092   return BRepAdaptor_Surface( TopoDS::Face( eos->_sWOL )).Value(_pos[i].X(), _pos[i].Y() ).XYZ();
8093 }
8094
8095 //================================================================================
8096 /*!
8097  * \brief Returns size and direction of the last segment
8098  */
8099 //================================================================================
8100
8101 gp_Ax1 _LayerEdge::LastSegment(double& segLen, _EdgesOnShape& eos) const
8102 {
8103   // find two non-coincident positions
8104   gp_XYZ orig = _pos.back();
8105   gp_XYZ vec;
8106   int iPrev = _pos.size() - 2;
8107   //const double tol = ( _len > 0 ) ? 0.3*_len : 1e-100; // adjusted for IPAL52478 + PAL22576
8108   const double tol = ( _len > 0 ) ? ( 1e-6 * _len ) : 1e-100;
8109   while ( iPrev >= 0 )
8110   {
8111     vec = orig - _pos[iPrev];
8112     if ( vec.SquareModulus() > tol*tol )
8113       break;
8114     else
8115       iPrev--;
8116   }
8117
8118   // make gp_Ax1
8119   gp_Ax1 segDir;
8120   if ( iPrev < 0 )
8121   {
8122     segDir.SetLocation( SMESH_TNodeXYZ( _nodes[0] ));
8123     segDir.SetDirection( _normal );
8124     segLen = 0;
8125   }
8126   else
8127   {
8128     gp_Pnt pPrev = _pos[ iPrev ];
8129     if ( !eos._sWOL.IsNull() )
8130     {
8131       TopLoc_Location loc;
8132       if ( eos.SWOLType() == TopAbs_EDGE )
8133       {
8134         double f,l;
8135         Handle(Geom_Curve) curve = BRep_Tool::Curve( TopoDS::Edge( eos._sWOL ), loc, f,l);
8136         pPrev = curve->Value( pPrev.X() ).Transformed( loc );
8137       }
8138       else
8139       {
8140         Handle(Geom_Surface) surface = BRep_Tool::Surface( TopoDS::Face( eos._sWOL ), loc );
8141         pPrev = surface->Value( pPrev.X(), pPrev.Y() ).Transformed( loc );
8142       }
8143       vec = SMESH_TNodeXYZ( _nodes.back() ) - pPrev.XYZ();
8144     }
8145     segDir.SetLocation( pPrev );
8146     segDir.SetDirection( vec );
8147     segLen = vec.Modulus();
8148   }
8149
8150   return segDir;
8151 }
8152
8153 //================================================================================
8154 /*!
8155  * \brief Return the last (or \a which) position of the target node on a FACE. 
8156  *  \param [in] F - the FACE this _LayerEdge is inflated along
8157  *  \param [in] which - index of position
8158  *  \return gp_XY - result UV
8159  */
8160 //================================================================================
8161
8162 gp_XY _LayerEdge::LastUV( const TopoDS_Face& F, _EdgesOnShape& eos, int which ) const
8163 {
8164   if ( F.IsSame( eos._sWOL )) // F is my FACE
8165     return gp_XY( _pos.back().X(), _pos.back().Y() );
8166
8167   if ( eos.SWOLType() != TopAbs_EDGE ) // wrong call
8168     return gp_XY( 1e100, 1e100 );
8169
8170   // _sWOL is EDGE of F; _pos.back().X() is the last U on the EDGE
8171   double f, l, u = _pos[ which < 0 ? _pos.size()-1 : which ].X();
8172   Handle(Geom2d_Curve) C2d = BRep_Tool::CurveOnSurface( TopoDS::Edge(eos._sWOL), F, f,l);
8173   if ( !C2d.IsNull() && f <= u && u <= l )
8174     return C2d->Value( u ).XY();
8175
8176   return gp_XY( 1e100, 1e100 );
8177 }
8178
8179 //================================================================================
8180 /*!
8181  * \brief Test intersection of the last segment with a given triangle
8182  *   using Moller-Trumbore algorithm
8183  * Intersection is detected if distance to intersection is less than _LayerEdge._len
8184  */
8185 //================================================================================
8186
8187 bool _LayerEdge::SegTriaInter( const gp_Ax1& lastSegment,
8188                                const gp_XYZ& vert0,
8189                                const gp_XYZ& vert1,
8190                                const gp_XYZ& vert2,
8191                                double&       t,
8192                                const double& EPSILON) const
8193 {
8194   const gp_Pnt& orig = lastSegment.Location();
8195   const gp_Dir& dir  = lastSegment.Direction();
8196
8197   /* calculate distance from vert0 to ray origin */
8198   //gp_XYZ tvec = orig.XYZ() - vert0;
8199
8200   //if ( tvec * dir > EPSILON )
8201     // intersected face is at back side of the temporary face this _LayerEdge belongs to
8202     //return false;
8203
8204   gp_XYZ edge1 = vert1 - vert0;
8205   gp_XYZ edge2 = vert2 - vert0;
8206
8207   /* begin calculating determinant - also used to calculate U parameter */
8208   gp_XYZ pvec = dir.XYZ() ^ edge2;
8209
8210   /* if determinant is near zero, ray lies in plane of triangle */
8211   double det = edge1 * pvec;
8212
8213   const double ANGL_EPSILON = 1e-12;
8214   if ( det > -ANGL_EPSILON && det < ANGL_EPSILON )
8215     return false;
8216
8217   /* calculate distance from vert0 to ray origin */
8218   gp_XYZ tvec = orig.XYZ() - vert0;
8219
8220   /* calculate U parameter and test bounds */
8221   double u = ( tvec * pvec ) / det;
8222   //if (u < 0.0 || u > 1.0)
8223   if ( u < -EPSILON || u > 1.0 + EPSILON )
8224     return false;
8225
8226   /* prepare to test V parameter */
8227   gp_XYZ qvec = tvec ^ edge1;
8228
8229   /* calculate V parameter and test bounds */
8230   double v = (dir.XYZ() * qvec) / det;
8231   //if ( v < 0.0 || u + v > 1.0 )
8232   if ( v < -EPSILON || u + v > 1.0 + EPSILON )
8233     return false;
8234
8235   /* calculate t, ray intersects triangle */
8236   t = (edge2 * qvec) / det;
8237
8238   //return true;
8239   return t > 0.;
8240 }
8241
8242 //================================================================================
8243 /*!
8244  * \brief _LayerEdge, located at a concave VERTEX of a FACE, moves target nodes of
8245  *        neighbor _LayerEdge's by it's own inflation vector.
8246  *  \param [in] eov - EOS of the VERTEX
8247  *  \param [in] eos - EOS of the FACE
8248  *  \param [in] step - inflation step
8249  *  \param [in,out] badSmooEdges - tangled _LayerEdge's
8250  */
8251 //================================================================================
8252
8253 void _LayerEdge::MoveNearConcaVer( const _EdgesOnShape*    eov,
8254                                    const _EdgesOnShape*    eos,
8255                                    const int               step,
8256                                    vector< _LayerEdge* > & badSmooEdges )
8257 {
8258   // check if any of _neibors is in badSmooEdges
8259   if ( std::find_first_of( _neibors.begin(), _neibors.end(),
8260                            badSmooEdges.begin(), badSmooEdges.end() ) == _neibors.end() )
8261     return;
8262
8263   // get all edges to move
8264
8265   set< _LayerEdge* > edges;
8266
8267   // find a distance between _LayerEdge on VERTEX and its neighbors
8268   gp_XYZ  curPosV = SMESH_TNodeXYZ( _nodes.back() );
8269   double dist2 = 0;
8270   for ( size_t i = 0; i < _neibors.size(); ++i )
8271   {
8272     _LayerEdge* nEdge = _neibors[i];
8273     if ( nEdge->_nodes[0]->getshapeId() == eos->_shapeID )
8274     {
8275       edges.insert( nEdge );
8276       dist2 = Max( dist2, ( curPosV - nEdge->_pos.back() ).SquareModulus() );
8277     }
8278   }
8279   // add _LayerEdge's close to curPosV
8280   size_t nbE;
8281   do {
8282     nbE = edges.size();
8283     for ( set< _LayerEdge* >::iterator e = edges.begin(); e != edges.end(); ++e )
8284     {
8285       _LayerEdge* edgeF = *e;
8286       for ( size_t i = 0; i < edgeF->_neibors.size(); ++i )
8287       {
8288         _LayerEdge* nEdge = edgeF->_neibors[i];
8289         if ( nEdge->_nodes[0]->getshapeId() == eos->_shapeID &&
8290              dist2 > ( curPosV - nEdge->_pos.back() ).SquareModulus() )
8291           edges.insert( nEdge );
8292       }
8293     }
8294   }
8295   while ( nbE < edges.size() );
8296
8297   // move the target node of the got edges
8298
8299   gp_XYZ prevPosV = PrevPos();
8300   if ( eov->SWOLType() == TopAbs_EDGE )
8301   {
8302     BRepAdaptor_Curve curve ( TopoDS::Edge( eov->_sWOL ));
8303     prevPosV = curve.Value( prevPosV.X() ).XYZ();
8304   }
8305   else if ( eov->SWOLType() == TopAbs_FACE )
8306   {
8307     BRepAdaptor_Surface surface( TopoDS::Face( eov->_sWOL ));
8308     prevPosV = surface.Value( prevPosV.X(), prevPosV.Y() ).XYZ();
8309   }
8310
8311   SMDS_FacePositionPtr fPos;
8312   //double r = 1. - Min( 0.9, step / 10. );
8313   for ( set< _LayerEdge* >::iterator e = edges.begin(); e != edges.end(); ++e )
8314   {
8315     _LayerEdge*       edgeF = *e;
8316     const gp_XYZ     prevVF = edgeF->PrevPos() - prevPosV;
8317     const gp_XYZ    newPosF = curPosV + prevVF;
8318     SMDS_MeshNode* tgtNodeF = const_cast<SMDS_MeshNode*>( edgeF->_nodes.back() );
8319     tgtNodeF->setXYZ( newPosF.X(), newPosF.Y(), newPosF.Z() );
8320     edgeF->_pos.back() = newPosF;
8321     dumpMoveComm( tgtNodeF, "MoveNearConcaVer" ); // debug
8322
8323     // set _curvature to make edgeF updated by putOnOffsetSurface()
8324     if ( !edgeF->_curvature )
8325       if (( fPos = edgeF->_nodes[0]->GetPosition() ))
8326       {
8327         edgeF->_curvature = new _Curvature;
8328         edgeF->_curvature->_r = 0;
8329         edgeF->_curvature->_k = 0;
8330         edgeF->_curvature->_h2lenRatio = 0;
8331         edgeF->_curvature->_uv.SetCoord( fPos->GetUParameter(), fPos->GetVParameter() );
8332       }
8333   }
8334   // gp_XYZ inflationVec( SMESH_TNodeXYZ( _nodes.back() ) -
8335   //                      SMESH_TNodeXYZ( _nodes[0]    ));
8336   // for ( set< _LayerEdge* >::iterator e = edges.begin(); e != edges.end(); ++e )
8337   // {
8338   //   _LayerEdge*      edgeF = *e;
8339   //   gp_XYZ          newPos = SMESH_TNodeXYZ( edgeF->_nodes[0] ) + inflationVec;
8340   //   SMDS_MeshNode* tgtNode = const_cast<SMDS_MeshNode*>( edgeF->_nodes.back() );
8341   //   tgtNode->setXYZ( newPos.X(), newPos.Y(), newPos.Z() );
8342   //   edgeF->_pos.back() = newPosF;
8343   //   dumpMoveComm( tgtNode, "MoveNearConcaVer" ); // debug
8344   // }
8345
8346   // smooth _LayerEdge's around moved nodes
8347   //size_t nbBadBefore = badSmooEdges.size();
8348   for ( set< _LayerEdge* >::iterator e = edges.begin(); e != edges.end(); ++e )
8349   {
8350     _LayerEdge* edgeF = *e;
8351     for ( size_t j = 0; j < edgeF->_neibors.size(); ++j )
8352       if ( edgeF->_neibors[j]->_nodes[0]->getshapeId() == eos->_shapeID )
8353         //&& !edges.count( edgeF->_neibors[j] ))
8354       {
8355         _LayerEdge* edgeFN = edgeF->_neibors[j];
8356         edgeFN->Unset( SMOOTHED );
8357         int nbBad = edgeFN->Smooth( step, /*isConcaFace=*/true, /*findBest=*/true );
8358         // if ( nbBad > 0 )
8359         // {
8360         //   gp_XYZ         newPos = SMESH_TNodeXYZ( edgeFN->_nodes[0] ) + inflationVec;
8361         //   const gp_XYZ& prevPos = edgeFN->_pos[ edgeFN->_pos.size()-2 ];
8362         //   int        nbBadAfter = edgeFN->_simplices.size();
8363         //   double vol;
8364         //   for ( size_t iS = 0; iS < edgeFN->_simplices.size(); ++iS )
8365         //   {
8366         //     nbBadAfter -= edgeFN->_simplices[iS].IsForward( &prevPos, &newPos, vol );
8367         //   }
8368         //   if ( nbBadAfter <= nbBad )
8369         //   {
8370         //     SMDS_MeshNode* tgtNode = const_cast<SMDS_MeshNode*>( edgeFN->_nodes.back() );
8371         //     tgtNode->setXYZ( newPos.X(), newPos.Y(), newPos.Z() );
8372         //     edgeF->_pos.back() = newPosF;
8373         //     dumpMoveComm( tgtNode, "MoveNearConcaVer 2" ); // debug
8374         //     nbBad = nbBadAfter;
8375         //   }
8376         // }
8377         if ( nbBad > 0 )
8378           badSmooEdges.push_back( edgeFN );
8379       }
8380   }
8381     // move a bit not smoothed around moved nodes
8382   //   for ( size_t i = nbBadBefore; i < badSmooEdges.size(); ++i )
8383   //   {
8384   //   _LayerEdge*      edgeF = badSmooEdges[i];
8385   //   SMDS_MeshNode* tgtNode = const_cast<SMDS_MeshNode*>( edgeF->_nodes.back() );
8386   //   gp_XYZ         newPos1 = SMESH_TNodeXYZ( edgeF->_nodes[0] ) + inflationVec;
8387   //   gp_XYZ         newPos2 = 0.5 * ( newPos1 + SMESH_TNodeXYZ( tgtNode ));
8388   //   tgtNode->setXYZ( newPos2.X(), newPos2.Y(), newPos2.Z() );
8389   //   edgeF->_pos.back() = newPosF;
8390   //   dumpMoveComm( tgtNode, "MoveNearConcaVer 2" ); // debug
8391   // }
8392 }
8393
8394 //================================================================================
8395 /*!
8396  * \brief Perform smooth of _LayerEdge's based on EDGE's
8397  *  \retval bool - true if node has been moved
8398  */
8399 //================================================================================
8400
8401 bool _LayerEdge::SmoothOnEdge(Handle(ShapeAnalysis_Surface)& surface,
8402                               const TopoDS_Face&             F,
8403                               SMESH_MesherHelper&            helper)
8404 {
8405   ASSERT( IsOnEdge() );
8406
8407   SMDS_MeshNode* tgtNode = const_cast<SMDS_MeshNode*>( _nodes.back() );
8408   SMESH_TNodeXYZ oldPos( tgtNode );
8409   double dist01, distNewOld;
8410   
8411   SMESH_TNodeXYZ p0( _2neibors->tgtNode(0));
8412   SMESH_TNodeXYZ p1( _2neibors->tgtNode(1));
8413   dist01 = p0.Distance( _2neibors->tgtNode(1) );
8414
8415   gp_Pnt newPos = p0 * _2neibors->_wgt[0] + p1 * _2neibors->_wgt[1];
8416   double lenDelta = 0;
8417   if ( _curvature )
8418   {
8419     //lenDelta = _curvature->lenDelta( _len );
8420     lenDelta = _curvature->lenDeltaByDist( dist01 );
8421     newPos.ChangeCoord() += _normal * lenDelta;
8422   }
8423
8424   distNewOld = newPos.Distance( oldPos );
8425
8426   if ( F.IsNull() )
8427   {
8428     if ( _2neibors->_plnNorm )
8429     {
8430       // put newPos on the plane defined by source node and _plnNorm
8431       gp_XYZ new2src = SMESH_TNodeXYZ( _nodes[0] ) - newPos.XYZ();
8432       double new2srcProj = (*_2neibors->_plnNorm) * new2src;
8433       newPos.ChangeCoord() += (*_2neibors->_plnNorm) * new2srcProj;
8434     }
8435     tgtNode->setXYZ( newPos.X(), newPos.Y(), newPos.Z() );
8436     _pos.back() = newPos.XYZ();
8437   }
8438   else
8439   {
8440     tgtNode->setXYZ( newPos.X(), newPos.Y(), newPos.Z() );
8441     gp_XY uv( Precision::Infinite(), 0 );
8442     helper.CheckNodeUV( F, tgtNode, uv, 1e-10, /*force=*/true );
8443     _pos.back().SetCoord( uv.X(), uv.Y(), 0 );
8444
8445     newPos = surface->Value( uv );
8446     tgtNode->setXYZ( newPos.X(), newPos.Y(), newPos.Z() );
8447   }
8448
8449   // commented for IPAL0052478
8450   // if ( _curvature && lenDelta < 0 )
8451   // {
8452   //   gp_Pnt prevPos( _pos[ _pos.size()-2 ]);
8453   //   _len -= prevPos.Distance( oldPos );
8454   //   _len += prevPos.Distance( newPos );
8455   // }
8456   bool moved = distNewOld > dist01/50;
8457   //if ( moved )
8458   dumpMove( tgtNode ); // debug
8459
8460   return moved;
8461 }
8462
8463 //================================================================================
8464 /*!
8465  * \brief Perform 3D smooth of nodes inflated from FACE. No check of validity
8466  */
8467 //================================================================================
8468
8469 void _LayerEdge::SmoothWoCheck()
8470 {
8471   if ( Is( DIFFICULT ))
8472     return;
8473
8474   bool moved = Is( SMOOTHED );
8475   for ( size_t i = 0; i < _neibors.size()  &&  !moved; ++i )
8476     moved = _neibors[i]->Is( SMOOTHED );
8477   if ( !moved )
8478     return;
8479
8480   gp_XYZ newPos = (this->*_smooFunction)(); // fun chosen by ChooseSmooFunction()
8481
8482   SMDS_MeshNode* n = const_cast< SMDS_MeshNode* >( _nodes.back() );
8483   n->setXYZ( newPos.X(), newPos.Y(), newPos.Z());
8484   _pos.back() = newPos;
8485
8486   dumpMoveComm( n, SMESH_Comment("No check - ") << _funNames[ smooFunID() ]);
8487 }
8488
8489 //================================================================================
8490 /*!
8491  * \brief Checks validity of _neibors on EDGEs and VERTEXes
8492  */
8493 //================================================================================
8494
8495 int _LayerEdge::CheckNeiborsOnBoundary( vector< _LayerEdge* >* badNeibors, bool * needSmooth )
8496 {
8497   if ( ! Is( NEAR_BOUNDARY ))
8498     return 0;
8499
8500   int nbBad = 0;
8501   double vol;
8502   for ( size_t iN = 0; iN < _neibors.size(); ++iN )
8503   {
8504     _LayerEdge* eN = _neibors[iN];
8505     if ( eN->_nodes[0]->getshapeId() == _nodes[0]->getshapeId() )
8506       continue;
8507     if ( needSmooth )
8508       *needSmooth |= ( eN->Is( _LayerEdge::BLOCKED ) ||
8509                        eN->Is( _LayerEdge::NORMAL_UPDATED ) ||
8510                        eN->_pos.size() != _pos.size() );
8511
8512     SMESH_TNodeXYZ curPosN ( eN->_nodes.back() );
8513     SMESH_TNodeXYZ prevPosN( eN->_nodes[0] );
8514     for ( size_t i = 0; i < eN->_simplices.size(); ++i )
8515       if ( eN->_nodes.size() > 1 &&
8516            eN->_simplices[i].Includes( _nodes.back() ) &&
8517            !eN->_simplices[i].IsForward( &prevPosN, &curPosN, vol ))
8518       {
8519         ++nbBad;
8520         if ( badNeibors )
8521         {
8522           badNeibors->push_back( eN );
8523           debugMsg("Bad boundary simplex ( "
8524                    << " "<< eN->_nodes[0]->GetID()
8525                    << " "<< eN->_nodes.back()->GetID()
8526                    << " "<< eN->_simplices[i]._nPrev->GetID()
8527                    << " "<< eN->_simplices[i]._nNext->GetID() << " )" );
8528         }
8529         else
8530         {
8531           break;
8532         }
8533       }
8534   }
8535   return nbBad;
8536 }
8537
8538 //================================================================================
8539 /*!
8540  * \brief Perform 'smart' 3D smooth of nodes inflated from FACE
8541  *  \retval int - nb of bad simplices around this _LayerEdge
8542  */
8543 //================================================================================
8544
8545 int _LayerEdge::Smooth(const int step, bool findBest, vector< _LayerEdge* >& toSmooth )
8546 {
8547   if ( !Is( MOVED ) || Is( SMOOTHED ) || Is( BLOCKED ))
8548     return 0; // shape of simplices not changed
8549   if ( _simplices.size() < 2 )
8550     return 0; // _LayerEdge inflated along EDGE or FACE
8551
8552   if ( Is( DIFFICULT )) // || Is( ON_CONCAVE_FACE )
8553     findBest = true;
8554
8555   const gp_XYZ& curPos  = _pos.back();
8556   const gp_XYZ& prevPos = _pos[0]; //PrevPos();
8557
8558   // quality metrics (orientation) of tetras around _tgtNode
8559   int nbOkBefore = 0;
8560   double vol, minVolBefore = 1e100;
8561   for ( size_t i = 0; i < _simplices.size(); ++i )
8562   {
8563     nbOkBefore += _simplices[i].IsForward( &prevPos, &curPos, vol );
8564     minVolBefore = Min( minVolBefore, vol );
8565   }
8566   int nbBad = _simplices.size() - nbOkBefore;
8567
8568   bool bndNeedSmooth = false;
8569   if ( nbBad == 0 )
8570     nbBad = CheckNeiborsOnBoundary( 0, & bndNeedSmooth );
8571   if ( nbBad > 0 )
8572     Set( DISTORTED );
8573
8574   // evaluate min angle
8575   if ( nbBad == 0 && !findBest && !bndNeedSmooth )
8576   {
8577     size_t nbGoodAngles = _simplices.size();
8578     double angle;
8579     for ( size_t i = 0; i < _simplices.size(); ++i )
8580     {
8581       if ( !_simplices[i].IsMinAngleOK( curPos, angle ) && angle > _minAngle )
8582         --nbGoodAngles;
8583     }
8584     if ( nbGoodAngles == _simplices.size() )
8585     {
8586       Unset( MOVED );
8587       return 0;
8588     }
8589   }
8590   if ( Is( ON_CONCAVE_FACE ))
8591     findBest = true;
8592
8593   if ( step % 2 == 0 )
8594     findBest = false;
8595
8596   if ( Is( ON_CONCAVE_FACE ) && !findBest ) // alternate FUN_CENTROIDAL and FUN_LAPLACIAN
8597   {
8598     if ( _smooFunction == _funs[ FUN_LAPLACIAN ] )
8599       _smooFunction = _funs[ FUN_CENTROIDAL ];
8600     else
8601       _smooFunction = _funs[ FUN_LAPLACIAN ];
8602   }
8603
8604   // compute new position for the last _pos using different _funs
8605   gp_XYZ newPos;
8606   bool moved = false;
8607   for ( int iFun = -1; iFun < theNbSmooFuns; ++iFun )
8608   {
8609     if ( iFun < 0 )
8610       newPos = (this->*_smooFunction)(); // fun chosen by ChooseSmooFunction()
8611     else if ( _funs[ iFun ] == _smooFunction )
8612       continue; // _smooFunction again
8613     else if ( step > 1 )
8614       newPos = (this->*_funs[ iFun ])(); // try other smoothing fun
8615     else
8616       break; // let "easy" functions improve elements around distorted ones
8617
8618     if ( _curvature )
8619     {
8620       double delta  = _curvature->lenDelta( _len );
8621       if ( delta > 0 )
8622         newPos += _normal * delta;
8623       else
8624       {
8625         double segLen = _normal * ( newPos - prevPos );
8626         if ( segLen + delta > 0 )
8627           newPos += _normal * delta;
8628       }
8629       // double segLenChange = _normal * ( curPos - newPos );
8630       // newPos += 0.5 * _normal * segLenChange;
8631     }
8632
8633     int nbOkAfter = 0;
8634     double minVolAfter = 1e100;
8635     for ( size_t i = 0; i < _simplices.size(); ++i )
8636     {
8637       nbOkAfter += _simplices[i].IsForward( &prevPos, &newPos, vol );
8638       minVolAfter = Min( minVolAfter, vol );
8639     }
8640     // get worse?
8641     if ( nbOkAfter < nbOkBefore )
8642       continue;
8643
8644     if (( findBest ) &&
8645         ( nbOkAfter == nbOkBefore ) &&
8646         ( minVolAfter <= minVolBefore ))
8647       continue;
8648
8649     nbBad        = _simplices.size() - nbOkAfter;
8650     minVolBefore = minVolAfter;
8651     nbOkBefore   = nbOkAfter;
8652     moved        = true;
8653
8654     SMDS_MeshNode* n = const_cast< SMDS_MeshNode* >( _nodes.back() );
8655     n->setXYZ( newPos.X(), newPos.Y(), newPos.Z());
8656     _pos.back() = newPos;
8657
8658     dumpMoveComm( n, SMESH_Comment( _funNames[ iFun < 0 ? smooFunID() : iFun ] )
8659                   << (nbBad ? " --BAD" : ""));
8660
8661     if ( iFun > -1 )
8662     {
8663       continue; // look for a better function
8664     }
8665
8666     if ( !findBest )
8667       break;
8668
8669   } // loop on smoothing functions
8670
8671   if ( moved ) // notify _neibors
8672   {
8673     Set( SMOOTHED );
8674     for ( size_t i = 0; i < _neibors.size(); ++i )
8675       if ( !_neibors[i]->Is( MOVED ))
8676       {
8677         _neibors[i]->Set( MOVED );
8678         toSmooth.push_back( _neibors[i] );
8679       }
8680   }
8681
8682   return nbBad;
8683 }
8684
8685 //================================================================================
8686 /*!
8687  * \brief Perform 'smart' 3D smooth of nodes inflated from FACE
8688  *  \retval int - nb of bad simplices around this _LayerEdge
8689  */
8690 //================================================================================
8691
8692 int _LayerEdge::Smooth(const int step, const bool isConcaveFace, bool findBest )
8693 {
8694   if ( !_smooFunction )
8695     return 0; // _LayerEdge inflated along EDGE or FACE
8696   if ( Is( BLOCKED ))
8697     return 0; // not inflated
8698
8699   const gp_XYZ& curPos  = _pos.back();
8700   const gp_XYZ& prevPos = _pos[0]; //PrevCheckPos();
8701
8702   // quality metrics (orientation) of tetras around _tgtNode
8703   int nbOkBefore = 0;
8704   double vol, minVolBefore = 1e100;
8705   for ( size_t i = 0; i < _simplices.size(); ++i )
8706   {
8707     nbOkBefore += _simplices[i].IsForward( &prevPos, &curPos, vol );
8708     minVolBefore = Min( minVolBefore, vol );
8709   }
8710   int nbBad = _simplices.size() - nbOkBefore;
8711
8712   if ( isConcaveFace ) // alternate FUN_CENTROIDAL and FUN_LAPLACIAN
8713   {
8714     if      ( _smooFunction == _funs[ FUN_CENTROIDAL ] && step % 2 )
8715       _smooFunction = _funs[ FUN_LAPLACIAN ];
8716     else if ( _smooFunction == _funs[ FUN_LAPLACIAN ] && !( step % 2 ))
8717       _smooFunction = _funs[ FUN_CENTROIDAL ];
8718   }
8719
8720   // compute new position for the last _pos using different _funs
8721   gp_XYZ newPos;
8722   for ( int iFun = -1; iFun < theNbSmooFuns; ++iFun )
8723   {
8724     if ( iFun < 0 )
8725       newPos = (this->*_smooFunction)(); // fun chosen by ChooseSmooFunction()
8726     else if ( _funs[ iFun ] == _smooFunction )
8727       continue; // _smooFunction again
8728     else if ( step > 1 )
8729       newPos = (this->*_funs[ iFun ])(); // try other smoothing fun
8730     else
8731       break; // let "easy" functions improve elements around distorted ones
8732
8733     if ( _curvature )
8734     {
8735       double delta  = _curvature->lenDelta( _len );
8736       if ( delta > 0 )
8737         newPos += _normal * delta;
8738       else
8739       {
8740         double segLen = _normal * ( newPos - prevPos );
8741         if ( segLen + delta > 0 )
8742           newPos += _normal * delta;
8743       }
8744       // double segLenChange = _normal * ( curPos - newPos );
8745       // newPos += 0.5 * _normal * segLenChange;
8746     }
8747
8748     int nbOkAfter = 0;
8749     double minVolAfter = 1e100;
8750     for ( size_t i = 0; i < _simplices.size(); ++i )
8751     {
8752       nbOkAfter += _simplices[i].IsForward( &prevPos, &newPos, vol );
8753       minVolAfter = Min( minVolAfter, vol );
8754     }
8755     // get worse?
8756     if ( nbOkAfter < nbOkBefore )
8757       continue;
8758     if (( isConcaveFace || findBest ) &&
8759         ( nbOkAfter == nbOkBefore ) &&
8760         ( minVolAfter <= minVolBefore )
8761         )
8762       continue;
8763
8764     nbBad        = _simplices.size() - nbOkAfter;
8765     minVolBefore = minVolAfter;
8766     nbOkBefore   = nbOkAfter;
8767
8768     SMDS_MeshNode* n = const_cast< SMDS_MeshNode* >( _nodes.back() );
8769     n->setXYZ( newPos.X(), newPos.Y(), newPos.Z());
8770     _pos.back() = newPos;
8771
8772     dumpMoveComm( n, SMESH_Comment( _funNames[ iFun < 0 ? smooFunID() : iFun ] )
8773                   << ( nbBad ? "--BAD" : ""));
8774
8775     // commented for IPAL0052478
8776     // _len -= prevPos.Distance(SMESH_TNodeXYZ( n ));
8777     // _len += prevPos.Distance(newPos);
8778
8779     if ( iFun > -1 ) // findBest || the chosen _fun makes worse
8780     {
8781       //_smooFunction = _funs[ iFun ];
8782       // cout << "# " << _funNames[ iFun ] << "\t N:" << _nodes.back()->GetID()
8783       // << "\t nbBad: " << _simplices.size() - nbOkAfter
8784       // << " minVol: " << minVolAfter
8785       // << " " << newPos.X() << " " << newPos.Y() << " " << newPos.Z()
8786       // << endl;
8787       continue; // look for a better function
8788     }
8789
8790     if ( !findBest )
8791       break;
8792
8793   } // loop on smoothing functions
8794
8795   return nbBad;
8796 }
8797
8798 //================================================================================
8799 /*!
8800  * \brief Chooses a smoothing technique giving a position most close to an initial one.
8801  *        For a correct result, _simplices must contain nodes lying on geometry.
8802  */
8803 //================================================================================
8804
8805 void _LayerEdge::ChooseSmooFunction( const set< TGeomID >& concaveVertices,
8806                                      const TNode2Edge&     n2eMap)
8807 {
8808   if ( _smooFunction ) return;
8809
8810   // use smoothNefPolygon() near concaveVertices
8811   if ( !concaveVertices.empty() )
8812   {
8813     _smooFunction = _funs[ FUN_CENTROIDAL ];
8814
8815     Set( ON_CONCAVE_FACE );
8816
8817     for ( size_t i = 0; i < _simplices.size(); ++i )
8818     {
8819       if ( concaveVertices.count( _simplices[i]._nPrev->getshapeId() ))
8820       {
8821         _smooFunction = _funs[ FUN_NEFPOLY ];
8822
8823         // set FUN_CENTROIDAL to neighbor edges
8824         for ( i = 0; i < _neibors.size(); ++i )
8825         {
8826           if ( _neibors[i]->_nodes[0]->GetPosition()->GetDim() == 2 )
8827           {
8828             _neibors[i]->_smooFunction = _funs[ FUN_CENTROIDAL ];
8829           }
8830         }
8831         return;
8832       }
8833     }
8834
8835     // // this choice is done only if ( !concaveVertices.empty() ) for Grids/smesh/bugs_19/X1
8836     // // where the nodes are smoothed too far along a sphere thus creating
8837     // // inverted _simplices
8838     // double dist[theNbSmooFuns];
8839     // //double coef[theNbSmooFuns] = { 1., 1.2, 1.4, 1.4 };
8840     // double coef[theNbSmooFuns] = { 1., 1., 1., 1. };
8841
8842     // double minDist = Precision::Infinite();
8843     // gp_Pnt p = SMESH_TNodeXYZ( _nodes[0] );
8844     // for ( int i = 0; i < FUN_NEFPOLY; ++i )
8845     // {
8846     //   gp_Pnt newP = (this->*_funs[i])();
8847     //   dist[i] = p.SquareDistance( newP );
8848     //   if ( dist[i]*coef[i] < minDist )
8849     //   {
8850     //     _smooFunction = _funs[i];
8851     //     minDist = dist[i]*coef[i];
8852     //   }
8853     // }
8854   }
8855   else
8856   {
8857     _smooFunction = _funs[ FUN_LAPLACIAN ];
8858   }
8859   // int minDim = 3;
8860   // for ( size_t i = 0; i < _simplices.size(); ++i )
8861   //   minDim = Min( minDim, _simplices[i]._nPrev->GetPosition()->GetDim() );
8862   // if ( minDim == 0 )
8863   //   _smooFunction = _funs[ FUN_CENTROIDAL ];
8864   // else if ( minDim == 1 )
8865   //   _smooFunction = _funs[ FUN_CENTROIDAL ];
8866
8867
8868   // int iMin;
8869   // for ( int i = 0; i < FUN_NB; ++i )
8870   // {
8871   //   //cout << dist[i] << " ";
8872   //   if ( _smooFunction == _funs[i] ) {
8873   //     iMin = i;
8874   //     //debugMsg( fNames[i] );
8875   //     break;
8876   //   }
8877   // }
8878   // cout << _funNames[ iMin ] << "\t N:" << _nodes.back()->GetID() << endl;
8879 }
8880
8881 //================================================================================
8882 /*!
8883  * \brief Returns a name of _SmooFunction
8884  */
8885 //================================================================================
8886
8887 int _LayerEdge::smooFunID( _LayerEdge::PSmooFun fun) const
8888 {
8889   if ( !fun )
8890     fun = _smooFunction;
8891   for ( int i = 0; i < theNbSmooFuns; ++i )
8892     if ( fun == _funs[i] )
8893       return i;
8894
8895   return theNbSmooFuns;
8896 }
8897
8898 //================================================================================
8899 /*!
8900  * \brief Computes a new node position using Laplacian smoothing
8901  */
8902 //================================================================================
8903
8904 gp_XYZ _LayerEdge::smoothLaplacian()
8905 {
8906   gp_XYZ newPos (0,0,0);
8907   for ( size_t i = 0; i < _simplices.size(); ++i )
8908     newPos += SMESH_TNodeXYZ( _simplices[i]._nPrev );
8909   newPos /= _simplices.size();
8910
8911   return newPos;
8912 }
8913
8914 //================================================================================
8915 /*!
8916  * \brief Computes a new node position using angular-based smoothing
8917  */
8918 //================================================================================
8919
8920 gp_XYZ _LayerEdge::smoothAngular()
8921 {
8922   vector< gp_Vec > edgeDir;  edgeDir. reserve( _simplices.size() + 1 );
8923   vector< double > edgeSize; edgeSize.reserve( _simplices.size()     );
8924   vector< gp_XYZ > points;   points.  reserve( _simplices.size() + 1 );
8925
8926   gp_XYZ pPrev = SMESH_TNodeXYZ( _simplices.back()._nPrev );
8927   gp_XYZ pN( 0,0,0 );
8928   for ( size_t i = 0; i < _simplices.size(); ++i )
8929   {
8930     gp_XYZ p = SMESH_TNodeXYZ( _simplices[i]._nPrev );
8931     edgeDir.push_back( p - pPrev );
8932     edgeSize.push_back( edgeDir.back().Magnitude() );
8933     if ( edgeSize.back() < numeric_limits<double>::min() )
8934     {
8935       edgeDir.pop_back();
8936       edgeSize.pop_back();
8937     }
8938     else
8939     {
8940       edgeDir.back() /= edgeSize.back();
8941       points.push_back( p );
8942       pN += p;
8943     }
8944     pPrev = p;
8945   }
8946   edgeDir.push_back ( edgeDir[0] );
8947   edgeSize.push_back( edgeSize[0] );
8948   pN /= points.size();
8949
8950   gp_XYZ newPos(0,0,0);
8951   double sumSize = 0;
8952   for ( size_t i = 0; i < points.size(); ++i )
8953   {
8954     gp_Vec toN    = pN - points[i];
8955     double toNLen = toN.Magnitude();
8956     if ( toNLen < numeric_limits<double>::min() )
8957     {
8958       newPos += pN;
8959       continue;
8960     }
8961     gp_Vec bisec    = edgeDir[i] + edgeDir[i+1];
8962     double bisecLen = bisec.SquareMagnitude();
8963     if ( bisecLen < numeric_limits<double>::min() )
8964     {
8965       gp_Vec norm = edgeDir[i] ^ toN;
8966       bisec = norm ^ edgeDir[i];
8967       bisecLen = bisec.SquareMagnitude();
8968     }
8969     bisecLen = Sqrt( bisecLen );
8970     bisec /= bisecLen;
8971
8972 #if 1
8973     gp_XYZ pNew = ( points[i] + bisec.XYZ() * toNLen ) * bisecLen;
8974     sumSize += bisecLen;
8975 #else
8976     gp_XYZ pNew = ( points[i] + bisec.XYZ() * toNLen ) * ( edgeSize[i] + edgeSize[i+1] );
8977     sumSize += ( edgeSize[i] + edgeSize[i+1] );
8978 #endif
8979     newPos += pNew;
8980   }
8981   newPos /= sumSize;
8982
8983   // project newPos to an average plane
8984
8985   gp_XYZ norm(0,0,0); // plane normal
8986   points.push_back( points[0] );
8987   for ( size_t i = 1; i < points.size(); ++i )
8988   {
8989     gp_XYZ vec1 = points[ i-1 ] - pN;
8990     gp_XYZ vec2 = points[ i   ] - pN;
8991     gp_XYZ cross = vec1 ^ vec2;
8992     try {
8993       cross.Normalize();
8994       if ( cross * norm < numeric_limits<double>::min() )
8995         norm += cross.Reversed();
8996       else
8997         norm += cross;
8998     }
8999     catch (Standard_Failure) { // if |cross| == 0.
9000     }
9001   }
9002   gp_XYZ vec = newPos - pN;
9003   double r   = ( norm * vec ) / norm.SquareModulus();  // param [0,1] on norm
9004   newPos     = newPos - r * norm;
9005
9006   return newPos;
9007 }
9008
9009 //================================================================================
9010 /*!
9011  * \brief Computes a new node position using weighted node positions
9012  */
9013 //================================================================================
9014
9015 gp_XYZ _LayerEdge::smoothLengthWeighted()
9016 {
9017   vector< double > edgeSize; edgeSize.reserve( _simplices.size() + 1);
9018   vector< gp_XYZ > points;   points.  reserve( _simplices.size() );
9019
9020   gp_XYZ pPrev = SMESH_TNodeXYZ( _simplices.back()._nPrev );
9021   for ( size_t i = 0; i < _simplices.size(); ++i )
9022   {
9023     gp_XYZ p = SMESH_TNodeXYZ( _simplices[i]._nPrev );
9024     edgeSize.push_back( ( p - pPrev ).Modulus() );
9025     if ( edgeSize.back() < numeric_limits<double>::min() )
9026     {
9027       edgeSize.pop_back();
9028     }
9029     else
9030     {
9031       points.push_back( p );
9032     }
9033     pPrev = p;
9034   }
9035   edgeSize.push_back( edgeSize[0] );
9036
9037   gp_XYZ newPos(0,0,0);
9038   double sumSize = 0;
9039   for ( size_t i = 0; i < points.size(); ++i )
9040   {
9041     newPos += points[i] * ( edgeSize[i] + edgeSize[i+1] );
9042     sumSize += edgeSize[i] + edgeSize[i+1];
9043   }
9044   newPos /= sumSize;
9045   return newPos;
9046 }
9047
9048 //================================================================================
9049 /*!
9050  * \brief Computes a new node position using angular-based smoothing
9051  */
9052 //================================================================================
9053
9054 gp_XYZ _LayerEdge::smoothCentroidal()
9055 {
9056   gp_XYZ newPos(0,0,0);
9057   gp_XYZ pN = SMESH_TNodeXYZ( _nodes.back() );
9058   double sumSize = 0;
9059   for ( size_t i = 0; i < _simplices.size(); ++i )
9060   {
9061     gp_XYZ p1 = SMESH_TNodeXYZ( _simplices[i]._nPrev );
9062     gp_XYZ p2 = SMESH_TNodeXYZ( _simplices[i]._nNext );
9063     gp_XYZ gc = ( pN + p1 + p2 ) / 3.;
9064     double size = (( p1 - pN ) ^ ( p2 - pN )).Modulus();
9065
9066     sumSize += size;
9067     newPos += gc * size;
9068   }
9069   newPos /= sumSize;
9070
9071   return newPos;
9072 }
9073
9074 //================================================================================
9075 /*!
9076  * \brief Computes a new node position located inside a Nef polygon
9077  */
9078 //================================================================================
9079
9080 gp_XYZ _LayerEdge::smoothNefPolygon()
9081 #ifdef OLD_NEF_POLYGON
9082 {
9083   gp_XYZ newPos(0,0,0);
9084
9085   // get a plane to search a solution on
9086
9087   vector< gp_XYZ > vecs( _simplices.size() + 1 );
9088   size_t i;
9089   const double tol = numeric_limits<double>::min();
9090   gp_XYZ center(0,0,0);
9091   for ( i = 0; i < _simplices.size(); ++i )
9092   {
9093     vecs[i] = ( SMESH_TNodeXYZ( _simplices[i]._nNext ) -
9094                 SMESH_TNodeXYZ( _simplices[i]._nPrev ));
9095     center += SMESH_TNodeXYZ( _simplices[i]._nPrev );
9096   }
9097   vecs.back() = vecs[0];
9098   center /= _simplices.size();
9099
9100   gp_XYZ zAxis(0,0,0);
9101   for ( i = 0; i < _simplices.size(); ++i )
9102     zAxis += vecs[i] ^ vecs[i+1];
9103
9104   gp_XYZ yAxis;
9105   for ( i = 0; i < _simplices.size(); ++i )
9106   {
9107     yAxis = vecs[i];
9108     if ( yAxis.SquareModulus() > tol )
9109       break;
9110   }
9111   gp_XYZ xAxis = yAxis ^ zAxis;
9112   // SMESH_TNodeXYZ p0( _simplices[0]._nPrev );
9113   // const double tol = 1e-6 * ( p0.Distance( _simplices[1]._nPrev ) +
9114   //                             p0.Distance( _simplices[2]._nPrev ));
9115   // gp_XYZ center = smoothLaplacian();
9116   // gp_XYZ xAxis, yAxis, zAxis;
9117   // for ( i = 0; i < _simplices.size(); ++i )
9118   // {
9119   //   xAxis = SMESH_TNodeXYZ( _simplices[i]._nPrev ) - center;
9120   //   if ( xAxis.SquareModulus() > tol*tol )
9121   //     break;
9122   // }
9123   // for ( i = 1; i < _simplices.size(); ++i )
9124   // {
9125   //   yAxis = SMESH_TNodeXYZ( _simplices[i]._nPrev ) - center;
9126   //   zAxis = xAxis ^ yAxis;
9127   //   if ( zAxis.SquareModulus() > tol*tol )
9128   //     break;
9129   // }
9130   // if ( i == _simplices.size() ) return newPos;
9131
9132   yAxis = zAxis ^ xAxis;
9133   xAxis /= xAxis.Modulus();
9134   yAxis /= yAxis.Modulus();
9135
9136   // get half-planes of _simplices
9137
9138   vector< _halfPlane > halfPlns( _simplices.size() );
9139   int nbHP = 0;
9140   for ( size_t i = 0; i < _simplices.size(); ++i )
9141   {
9142     gp_XYZ OP1 = SMESH_TNodeXYZ( _simplices[i]._nPrev ) - center;
9143     gp_XYZ OP2 = SMESH_TNodeXYZ( _simplices[i]._nNext ) - center;
9144     gp_XY  p1( OP1 * xAxis, OP1 * yAxis );
9145     gp_XY  p2( OP2 * xAxis, OP2 * yAxis );
9146     gp_XY  vec12 = p2 - p1;
9147     double dist12 = vec12.Modulus();
9148     if ( dist12 < tol )
9149       continue;
9150     vec12 /= dist12;
9151     halfPlns[ nbHP ]._pos = p1;
9152     halfPlns[ nbHP ]._dir = vec12;
9153     halfPlns[ nbHP ]._inNorm.SetCoord( -vec12.Y(), vec12.X() );
9154     ++nbHP;
9155   }
9156
9157   // intersect boundaries of half-planes, define state of intersection points
9158   // in relation to all half-planes and calculate internal point of a 2D polygon
9159
9160   double sumLen = 0;
9161   gp_XY newPos2D (0,0);
9162
9163   enum { UNDEF = -1, NOT_OUT, IS_OUT, NO_INT };
9164   typedef std::pair< gp_XY, int > TIntPntState; // coord and isOut state
9165   TIntPntState undefIPS( gp_XY(1e100,1e100), UNDEF );
9166
9167   vector< vector< TIntPntState > > allIntPnts( nbHP );
9168   for ( int iHP1 = 0; iHP1 < nbHP; ++iHP1 )
9169   {
9170     vector< TIntPntState > & intPnts1 = allIntPnts[ iHP1 ];
9171     if ( intPnts1.empty() ) intPnts1.resize( nbHP, undefIPS );
9172
9173     int iPrev = SMESH_MesherHelper::WrapIndex( iHP1 - 1, nbHP );
9174     int iNext = SMESH_MesherHelper::WrapIndex( iHP1 + 1, nbHP );
9175
9176     int nbNotOut = 0;
9177     const gp_XY* segEnds[2] = { 0, 0 }; // NOT_OUT points
9178
9179     for ( int iHP2 = 0; iHP2 < nbHP; ++iHP2 )
9180     {
9181       if ( iHP1 == iHP2 ) continue;
9182
9183       TIntPntState & ips1 = intPnts1[ iHP2 ];
9184       if ( ips1.second == UNDEF )
9185       {
9186         // find an intersection point of boundaries of iHP1 and iHP2
9187
9188         if ( iHP2 == iPrev ) // intersection with neighbors is known
9189           ips1.first = halfPlns[ iHP1 ]._pos;
9190         else if ( iHP2 == iNext )
9191           ips1.first = halfPlns[ iHP2 ]._pos;
9192         else if ( !halfPlns[ iHP1 ].FindIntersection( halfPlns[ iHP2 ], ips1.first ))
9193           ips1.second = NO_INT;
9194
9195         // classify the found intersection point
9196         if ( ips1.second != NO_INT )
9197         {
9198           ips1.second = NOT_OUT;
9199           for ( int i = 0; i < nbHP && ips1.second == NOT_OUT; ++i )
9200             if ( i != iHP1 && i != iHP2 &&
9201                  halfPlns[ i ].IsOut( ips1.first, tol ))
9202               ips1.second = IS_OUT;
9203         }
9204         vector< TIntPntState > & intPnts2 = allIntPnts[ iHP2 ];
9205         if ( intPnts2.empty() ) intPnts2.resize( nbHP, undefIPS );
9206         TIntPntState & ips2 = intPnts2[ iHP1 ];
9207         ips2 = ips1;
9208       }
9209       if ( ips1.second == NOT_OUT )
9210       {
9211         ++nbNotOut;
9212         segEnds[ bool(segEnds[0]) ] = & ips1.first;
9213       }
9214     }
9215
9216     // find a NOT_OUT segment of boundary which is located between
9217     // two NOT_OUT int points
9218
9219     if ( nbNotOut < 2 )
9220       continue; // no such a segment
9221
9222     if ( nbNotOut > 2 )
9223     {
9224       // sort points along the boundary
9225       map< double, TIntPntState* > ipsByParam;
9226       for ( int iHP2 = 0; iHP2 < nbHP; ++iHP2 )
9227       {
9228         TIntPntState & ips1 = intPnts1[ iHP2 ];
9229         if ( ips1.second != NO_INT )
9230         {
9231           gp_XY     op = ips1.first - halfPlns[ iHP1 ]._pos;
9232           double param = op * halfPlns[ iHP1 ]._dir;
9233           ipsByParam.insert( make_pair( param, & ips1 ));
9234         }
9235       }
9236       // look for two neighboring NOT_OUT points
9237       nbNotOut = 0;
9238       map< double, TIntPntState* >::iterator u2ips = ipsByParam.begin();
9239       for ( ; u2ips != ipsByParam.end(); ++u2ips )
9240       {
9241         TIntPntState & ips1 = *(u2ips->second);
9242         if ( ips1.second == NOT_OUT )
9243           segEnds[ bool( nbNotOut++ ) ] = & ips1.first;
9244         else if ( nbNotOut >= 2 )
9245           break;
9246         else
9247           nbNotOut = 0;
9248       }
9249     }
9250
9251     if ( nbNotOut >= 2 )
9252     {
9253       double len = ( *segEnds[0] - *segEnds[1] ).Modulus();
9254       sumLen += len;
9255
9256       newPos2D += 0.5 * len * ( *segEnds[0] + *segEnds[1] );
9257     }
9258   }
9259
9260   if ( sumLen > 0 )
9261   {
9262     newPos2D /= sumLen;
9263     newPos = center + xAxis * newPos2D.X() + yAxis * newPos2D.Y();
9264   }
9265   else
9266   {
9267     newPos = center;
9268   }
9269
9270   return newPos;
9271 }
9272 #else // OLD_NEF_POLYGON
9273 { ////////////////////////////////// NEW
9274   gp_XYZ newPos(0,0,0);
9275
9276   // get a plane to search a solution on
9277
9278   size_t i;
9279   gp_XYZ center(0,0,0);
9280   for ( i = 0; i < _simplices.size(); ++i )
9281     center += SMESH_TNodeXYZ( _simplices[i]._nPrev );
9282   center /= _simplices.size();
9283
9284   vector< gp_XYZ > vecs( _simplices.size() + 1 );
9285   for ( i = 0; i < _simplices.size(); ++i )
9286     vecs[i] = SMESH_TNodeXYZ( _simplices[i]._nPrev ) - center;
9287   vecs.back() = vecs[0];
9288
9289   const double tol = numeric_limits<double>::min();
9290   gp_XYZ zAxis(0,0,0);
9291   for ( i = 0; i < _simplices.size(); ++i )
9292   {
9293     gp_XYZ cross = vecs[i] ^ vecs[i+1];
9294     try {
9295       cross.Normalize();
9296       if ( cross * zAxis < tol )
9297         zAxis += cross.Reversed();
9298       else
9299         zAxis += cross;
9300     }
9301     catch (Standard_Failure) { // if |cross| == 0.
9302     }
9303   }
9304
9305   gp_XYZ yAxis;
9306   for ( i = 0; i < _simplices.size(); ++i )
9307   {
9308     yAxis = vecs[i];
9309     if ( yAxis.SquareModulus() > tol )
9310       break;
9311   }
9312   gp_XYZ xAxis = yAxis ^ zAxis;
9313   // SMESH_TNodeXYZ p0( _simplices[0]._nPrev );
9314   // const double tol = 1e-6 * ( p0.Distance( _simplices[1]._nPrev ) +
9315   //                             p0.Distance( _simplices[2]._nPrev ));
9316   // gp_XYZ center = smoothLaplacian();
9317   // gp_XYZ xAxis, yAxis, zAxis;
9318   // for ( i = 0; i < _simplices.size(); ++i )
9319   // {
9320   //   xAxis = SMESH_TNodeXYZ( _simplices[i]._nPrev ) - center;
9321   //   if ( xAxis.SquareModulus() > tol*tol )
9322   //     break;
9323   // }
9324   // for ( i = 1; i < _simplices.size(); ++i )
9325   // {
9326   //   yAxis = SMESH_TNodeXYZ( _simplices[i]._nPrev ) - center;
9327   //   zAxis = xAxis ^ yAxis;
9328   //   if ( zAxis.SquareModulus() > tol*tol )
9329   //     break;
9330   // }
9331   // if ( i == _simplices.size() ) return newPos;
9332
9333   yAxis = zAxis ^ xAxis;
9334   xAxis /= xAxis.Modulus();
9335   yAxis /= yAxis.Modulus();
9336
9337   // get half-planes of _simplices
9338
9339   vector< _halfPlane > halfPlns( _simplices.size() );
9340   int nbHP = 0;
9341   for ( size_t i = 0; i < _simplices.size(); ++i )
9342   {
9343     const gp_XYZ& OP1 = vecs[ i   ];
9344     const gp_XYZ& OP2 = vecs[ i+1 ];
9345     gp_XY  p1( OP1 * xAxis, OP1 * yAxis );
9346     gp_XY  p2( OP2 * xAxis, OP2 * yAxis );
9347     gp_XY  vec12 = p2 - p1;
9348     double dist12 = vec12.Modulus();
9349     if ( dist12 < tol )
9350       continue;
9351     vec12 /= dist12;
9352     halfPlns[ nbHP ]._pos = p1;
9353     halfPlns[ nbHP ]._dir = vec12;
9354     halfPlns[ nbHP ]._inNorm.SetCoord( -vec12.Y(), vec12.X() );
9355     ++nbHP;
9356   }
9357
9358   // intersect boundaries of half-planes, define state of intersection points
9359   // in relation to all half-planes and calculate internal point of a 2D polygon
9360
9361   double sumLen = 0;
9362   gp_XY newPos2D (0,0);
9363
9364   enum { UNDEF = -1, NOT_OUT, IS_OUT, NO_INT };
9365   typedef std::pair< gp_XY, int > TIntPntState; // coord and isOut state
9366   TIntPntState undefIPS( gp_XY(1e100,1e100), UNDEF );
9367
9368   vector< vector< TIntPntState > > allIntPnts( nbHP );
9369   for ( int iHP1 = 0; iHP1 < nbHP; ++iHP1 )
9370   {
9371     vector< TIntPntState > & intPnts1 = allIntPnts[ iHP1 ];
9372     if ( intPnts1.empty() ) intPnts1.resize( nbHP, undefIPS );
9373
9374     int iPrev = SMESH_MesherHelper::WrapIndex( iHP1 - 1, nbHP );
9375     int iNext = SMESH_MesherHelper::WrapIndex( iHP1 + 1, nbHP );
9376
9377     int nbNotOut = 0;
9378     const gp_XY* segEnds[2] = { 0, 0 }; // NOT_OUT points
9379
9380     for ( int iHP2 = 0; iHP2 < nbHP; ++iHP2 )
9381     {
9382       if ( iHP1 == iHP2 ) continue;
9383
9384       TIntPntState & ips1 = intPnts1[ iHP2 ];
9385       if ( ips1.second == UNDEF )
9386       {
9387         // find an intersection point of boundaries of iHP1 and iHP2
9388
9389         if ( iHP2 == iPrev ) // intersection with neighbors is known
9390           ips1.first = halfPlns[ iHP1 ]._pos;
9391         else if ( iHP2 == iNext )
9392           ips1.first = halfPlns[ iHP2 ]._pos;
9393         else if ( !halfPlns[ iHP1 ].FindIntersection( halfPlns[ iHP2 ], ips1.first ))
9394           ips1.second = NO_INT;
9395
9396         // classify the found intersection point
9397         if ( ips1.second != NO_INT )
9398         {
9399           ips1.second = NOT_OUT;
9400           for ( int i = 0; i < nbHP && ips1.second == NOT_OUT; ++i )
9401             if ( i != iHP1 && i != iHP2 &&
9402                  halfPlns[ i ].IsOut( ips1.first, tol ))
9403               ips1.second = IS_OUT;
9404         }
9405         vector< TIntPntState > & intPnts2 = allIntPnts[ iHP2 ];
9406         if ( intPnts2.empty() ) intPnts2.resize( nbHP, undefIPS );
9407         TIntPntState & ips2 = intPnts2[ iHP1 ];
9408         ips2 = ips1;
9409       }
9410       if ( ips1.second == NOT_OUT )
9411       {
9412         ++nbNotOut;
9413         segEnds[ bool(segEnds[0]) ] = & ips1.first;
9414       }
9415     }
9416
9417     // find a NOT_OUT segment of boundary which is located between
9418     // two NOT_OUT int points
9419
9420     if ( nbNotOut < 2 )
9421       continue; // no such a segment
9422
9423     if ( nbNotOut > 2 )
9424     {
9425       // sort points along the boundary
9426       map< double, TIntPntState* > ipsByParam;
9427       for ( int iHP2 = 0; iHP2 < nbHP; ++iHP2 )
9428       {
9429         TIntPntState & ips1 = intPnts1[ iHP2 ];
9430         if ( ips1.second != NO_INT )
9431         {
9432           gp_XY     op = ips1.first - halfPlns[ iHP1 ]._pos;
9433           double param = op * halfPlns[ iHP1 ]._dir;
9434           ipsByParam.insert( make_pair( param, & ips1 ));
9435         }
9436       }
9437       // look for two neighboring NOT_OUT points
9438       nbNotOut = 0;
9439       map< double, TIntPntState* >::iterator u2ips = ipsByParam.begin();
9440       for ( ; u2ips != ipsByParam.end(); ++u2ips )
9441       {
9442         TIntPntState & ips1 = *(u2ips->second);
9443         if ( ips1.second == NOT_OUT )
9444           segEnds[ bool( nbNotOut++ ) ] = & ips1.first;
9445         else if ( nbNotOut >= 2 )
9446           break;
9447         else
9448           nbNotOut = 0;
9449       }
9450     }
9451
9452     if ( nbNotOut >= 2 )
9453     {
9454       double len = ( *segEnds[0] - *segEnds[1] ).Modulus();
9455       sumLen += len;
9456
9457       newPos2D += 0.5 * len * ( *segEnds[0] + *segEnds[1] );
9458     }
9459   }
9460
9461   if ( sumLen > 0 )
9462   {
9463     newPos2D /= sumLen;
9464     newPos = center + xAxis * newPos2D.X() + yAxis * newPos2D.Y();
9465   }
9466   else
9467   {
9468     newPos = center;
9469   }
9470
9471   return newPos;
9472 }
9473 #endif // OLD_NEF_POLYGON
9474
9475 //================================================================================
9476 /*!
9477  * \brief Add a new segment to _LayerEdge during inflation
9478  */
9479 //================================================================================
9480
9481 void _LayerEdge::SetNewLength( double len, _EdgesOnShape& eos, SMESH_MesherHelper& helper )
9482 {
9483   if ( Is( BLOCKED ))
9484     return;
9485
9486   if ( len > _maxLen )
9487   {
9488     len = _maxLen;
9489     Block( eos.GetData() );
9490   }
9491   const double lenDelta = len - _len;
9492   if ( lenDelta < len * 1e-3  )
9493   {
9494     Block( eos.GetData() );
9495     return;
9496   }
9497
9498   SMDS_MeshNode* n = const_cast< SMDS_MeshNode*>( _nodes.back() );
9499   gp_XYZ oldXYZ = SMESH_TNodeXYZ( n );
9500   gp_XYZ newXYZ;
9501   if ( eos._hyp.IsOffsetMethod() )
9502   {
9503     newXYZ = oldXYZ;
9504     gp_Vec faceNorm;
9505     SMDS_ElemIteratorPtr faceIt = _nodes[0]->GetInverseElementIterator( SMDSAbs_Face );
9506     while ( faceIt->more() )
9507     {
9508       const SMDS_MeshElement* face = faceIt->next();
9509       if ( !eos.GetNormal( face, faceNorm ))
9510         continue;
9511
9512       // translate plane of a face
9513       gp_XYZ baryCenter = oldXYZ + faceNorm.XYZ() * lenDelta;
9514
9515       // find point of intersection of the face plane located at baryCenter
9516       // and _normal located at newXYZ
9517       double d   = -( faceNorm.XYZ() * baryCenter ); // d of plane equation ax+by+cz+d=0
9518       double dot =  ( faceNorm.XYZ() * _normal );
9519       if ( dot < std::numeric_limits<double>::min() )
9520         dot = lenDelta * 1e-3;
9521       double step = -( faceNorm.XYZ() * newXYZ + d ) / dot;
9522       newXYZ += step * _normal;
9523     }
9524     _lenFactor = _normal * ( newXYZ - oldXYZ ) / lenDelta; // _lenFactor is used in InvalidateStep()
9525   }
9526   else
9527   {
9528     newXYZ = oldXYZ + _normal * lenDelta * _lenFactor;
9529   }
9530
9531   n->setXYZ( newXYZ.X(), newXYZ.Y(), newXYZ.Z() );
9532   _pos.push_back( newXYZ );
9533
9534   if ( !eos._sWOL.IsNull() )
9535   {
9536     double distXYZ[4];
9537     bool uvOK = false;
9538     if ( eos.SWOLType() == TopAbs_EDGE )
9539     {
9540       double u = Precision::Infinite(); // to force projection w/o distance check
9541       uvOK = helper.CheckNodeU( TopoDS::Edge( eos._sWOL ), n, u,
9542                                 /*tol=*/2*lenDelta, /*force=*/true, distXYZ );
9543       _pos.back().SetCoord( u, 0, 0 );
9544       if ( _nodes.size() > 1 && uvOK )
9545       {
9546         SMDS_EdgePositionPtr pos = n->GetPosition();
9547         pos->SetUParameter( u );
9548       }
9549     }
9550     else //  TopAbs_FACE
9551     {
9552       gp_XY uv( Precision::Infinite(), 0 );
9553       uvOK = helper.CheckNodeUV( TopoDS::Face( eos._sWOL ), n, uv,
9554                                  /*tol=*/2*lenDelta, /*force=*/true, distXYZ );
9555       _pos.back().SetCoord( uv.X(), uv.Y(), 0 );
9556       if ( _nodes.size() > 1 && uvOK )
9557       {
9558         SMDS_FacePositionPtr pos = n->GetPosition();
9559         pos->SetUParameter( uv.X() );
9560         pos->SetVParameter( uv.Y() );
9561       }
9562     }
9563     if ( uvOK )
9564     {
9565       n->setXYZ( distXYZ[1], distXYZ[2], distXYZ[3]);
9566     }
9567     else
9568     {
9569       n->setXYZ( oldXYZ.X(), oldXYZ.Y(), oldXYZ.Z() );
9570       _pos.pop_back();
9571       Block( eos.GetData() );
9572       return;
9573     }
9574   }
9575
9576   _len = len;
9577
9578   // notify _neibors
9579   if ( eos.ShapeType() != TopAbs_FACE )
9580   {
9581     for ( size_t i = 0; i < _neibors.size(); ++i )
9582       //if (  _len > _neibors[i]->GetSmooLen() )
9583         _neibors[i]->Set( MOVED );
9584
9585     Set( MOVED );
9586   }
9587   dumpMove( n ); //debug
9588 }
9589
9590 //================================================================================
9591 /*!
9592  * \brief Set BLOCKED flag and propagate limited _maxLen to _neibors
9593  */
9594 //================================================================================
9595
9596 void _LayerEdge::Block( _SolidData& data )
9597 {
9598   //if ( Is( BLOCKED )) return;
9599   Set( BLOCKED );
9600
9601   SMESH_Comment msg( "#BLOCK shape=");
9602   msg << data.GetShapeEdges( this )->_shapeID
9603       << ", nodes " << _nodes[0]->GetID() << ", " << _nodes.back()->GetID();
9604   dumpCmd( msg + " -- BEGIN");
9605
9606   SetMaxLen( _len );
9607   std::queue<_LayerEdge*> queue;
9608   queue.push( this );
9609
9610   gp_Pnt pSrc, pTgt, pSrcN, pTgtN;
9611   while ( !queue.empty() )
9612   {
9613     _LayerEdge* edge = queue.front(); queue.pop();
9614     pSrc = SMESH_TNodeXYZ( edge->_nodes[0] );
9615     pTgt = SMESH_TNodeXYZ( edge->_nodes.back() );
9616     for ( size_t iN = 0; iN < edge->_neibors.size(); ++iN )
9617     {
9618       _LayerEdge* neibor = edge->_neibors[iN];
9619       if ( neibor->_maxLen < edge->_maxLen * 1.01 )
9620         continue;
9621       pSrcN = SMESH_TNodeXYZ( neibor->_nodes[0] );
9622       pTgtN = SMESH_TNodeXYZ( neibor->_nodes.back() );
9623       double minDist = pSrc.SquareDistance( pSrcN );
9624       minDist   = Min( pTgt.SquareDistance( pTgtN ), minDist );
9625       minDist   = Min( pSrc.SquareDistance( pTgtN ), minDist );
9626       minDist   = Min( pTgt.SquareDistance( pSrcN ), minDist );
9627       double newMaxLen = edge->_maxLen + 0.5 * Sqrt( minDist );
9628       //if ( edge->_nodes[0]->getshapeId() == neibor->_nodes[0]->getshapeId() ) viscous_layers_00/A3
9629       {
9630         //newMaxLen *= edge->_lenFactor / neibor->_lenFactor;
9631         // newMaxLen *= Min( edge->_lenFactor / neibor->_lenFactor,
9632         //                   neibor->_lenFactor / edge->_lenFactor );
9633       }
9634       if ( neibor->_maxLen > newMaxLen )
9635       {
9636         neibor->SetMaxLen( newMaxLen );
9637         if ( neibor->_maxLen < neibor->_len )
9638         {
9639           _EdgesOnShape* eos = data.GetShapeEdges( neibor );
9640           int       lastStep = neibor->Is( BLOCKED ) ? 1 : 0;
9641           while ( neibor->_len > neibor->_maxLen &&
9642                   neibor->NbSteps() > lastStep )
9643             neibor->InvalidateStep( neibor->NbSteps(), *eos, /*restoreLength=*/true );
9644           neibor->SetNewLength( neibor->_maxLen, *eos, data.GetHelper() );
9645           //neibor->Block( data );
9646         }
9647         queue.push( neibor );
9648       }
9649     }
9650   }
9651   dumpCmd( msg + " -- END");
9652 }
9653
9654 //================================================================================
9655 /*!
9656  * \brief Remove last inflation step
9657  */
9658 //================================================================================
9659
9660 void _LayerEdge::InvalidateStep( size_t curStep, const _EdgesOnShape& eos, bool restoreLength )
9661 {
9662   if ( _pos.size() > curStep && _nodes.size() > 1 )
9663   {
9664     _pos.resize( curStep );
9665
9666     gp_Pnt      nXYZ = _pos.back();
9667     SMDS_MeshNode* n = const_cast< SMDS_MeshNode*>( _nodes.back() );
9668     SMESH_TNodeXYZ curXYZ( n );
9669     if ( !eos._sWOL.IsNull() )
9670     {
9671       TopLoc_Location loc;
9672       if ( eos.SWOLType() == TopAbs_EDGE )
9673       {
9674         SMDS_EdgePositionPtr pos = n->GetPosition();
9675         pos->SetUParameter( nXYZ.X() );
9676         double f,l;
9677         Handle(Geom_Curve) curve = BRep_Tool::Curve( TopoDS::Edge( eos._sWOL ), loc, f,l);
9678         nXYZ = curve->Value( nXYZ.X() ).Transformed( loc );
9679       }
9680       else
9681       {
9682         SMDS_FacePositionPtr pos = n->GetPosition();
9683         pos->SetUParameter( nXYZ.X() );
9684         pos->SetVParameter( nXYZ.Y() );
9685         Handle(Geom_Surface) surface = BRep_Tool::Surface( TopoDS::Face(eos._sWOL), loc );
9686         nXYZ = surface->Value( nXYZ.X(), nXYZ.Y() ).Transformed( loc );
9687       }
9688     }
9689     n->setXYZ( nXYZ.X(), nXYZ.Y(), nXYZ.Z() );
9690     dumpMove( n );
9691
9692     if ( restoreLength )
9693     {
9694       if ( NbSteps() == 0 )
9695         _len = 0.;
9696       else if ( IsOnFace() && Is( MOVED ))
9697         _len = ( nXYZ.XYZ() - SMESH_NodeXYZ( _nodes[0] )) * _normal;
9698       else
9699         _len -= ( nXYZ.XYZ() - curXYZ ).Modulus() / _lenFactor;
9700     }
9701   }
9702   return;
9703 }
9704
9705 //================================================================================
9706 /*!
9707  * \brief Return index of a _pos distant from _normal
9708  */
9709 //================================================================================
9710
9711 int _LayerEdge::GetSmoothedPos( const double tol )
9712 {
9713   int iSmoothed = 0;
9714   for ( size_t i = 1; i < _pos.size() && !iSmoothed; ++i )
9715   {
9716     double normDist = ( _pos[i] - _pos[0] ).Crossed( _normal ).SquareModulus();
9717     if ( normDist > tol * tol )
9718       iSmoothed = i;
9719   }
9720   return iSmoothed;
9721 }
9722
9723 //================================================================================
9724 /*!
9725  * \brief Smooth a path formed by _pos of a _LayerEdge smoothed on FACE
9726  */
9727 //================================================================================
9728
9729 void _LayerEdge::SmoothPos( const vector< double >& segLen, const double tol )
9730 {
9731   if ( /*Is( NORMAL_UPDATED ) ||*/ _pos.size() <= 2 )
9732     return;
9733
9734   // find the 1st smoothed _pos
9735   int iSmoothed = GetSmoothedPos( tol );
9736   if ( !iSmoothed ) return;
9737
9738   gp_XYZ normal = _normal;
9739   if ( Is( NORMAL_UPDATED ))
9740   {
9741     double minDot = 1;
9742     for ( size_t i = 0; i < _neibors.size(); ++i )
9743     {
9744       if ( _neibors[i]->IsOnFace() )
9745       {
9746         double dot = _normal * _neibors[i]->_normal;
9747         if ( dot < minDot )
9748         {
9749           normal = _neibors[i]->_normal;
9750           minDot = dot;
9751         }
9752       }
9753     }
9754     if ( minDot == 1. )
9755       for ( size_t i = 1; i < _pos.size(); ++i )
9756       {
9757         normal = _pos[i] - _pos[0];
9758         double size = normal.Modulus();
9759         if ( size > RealSmall() )
9760         {
9761           normal /= size;
9762           break;
9763         }
9764       }
9765   }
9766   const double r = 0.2;
9767   for ( int iter = 0; iter < 50; ++iter )
9768   {
9769     double minDot = 1;
9770     for ( size_t i = Max( 1, iSmoothed-1-iter ); i < _pos.size()-1; ++i )
9771     {
9772       gp_XYZ midPos = 0.5 * ( _pos[i-1] + _pos[i+1] );
9773       gp_XYZ newPos = ( 1-r ) * midPos + r * _pos[i];
9774       _pos[i] = newPos;
9775       double midLen = 0.5 * ( segLen[i-1] + segLen[i+1] );
9776       double newLen = ( 1-r ) * midLen + r * segLen[i];
9777       const_cast< double& >( segLen[i] ) = newLen;
9778       // check angle between normal and (_pos[i+1], _pos[i] )
9779       gp_XYZ posDir = _pos[i+1] - _pos[i];
9780       double size   = posDir.SquareModulus();
9781       if ( size > RealSmall() )
9782         minDot = Min( minDot, ( normal * posDir ) * ( normal * posDir ) / size );
9783     }
9784     if ( minDot > 0.5 * 0.5 )
9785       break;
9786   }
9787   return;
9788 }
9789
9790 //================================================================================
9791 /*!
9792  * \brief Print flags
9793  */
9794 //================================================================================
9795
9796 std::string _LayerEdge::DumpFlags() const
9797 {
9798   SMESH_Comment dump;
9799   for ( int flag = 1; flag < 0x1000000; flag *= 2 )
9800     if ( _flags & flag )
9801     {
9802       EFlags f = (EFlags) flag;
9803       switch ( f ) {
9804       case TO_SMOOTH:       dump << "TO_SMOOTH";       break;
9805       case MOVED:           dump << "MOVED";           break;
9806       case SMOOTHED:        dump << "SMOOTHED";        break;
9807       case DIFFICULT:       dump << "DIFFICULT";       break;
9808       case ON_CONCAVE_FACE: dump << "ON_CONCAVE_FACE"; break;
9809       case BLOCKED:         dump << "BLOCKED";         break;
9810       case INTERSECTED:     dump << "INTERSECTED";     break;
9811       case NORMAL_UPDATED:  dump << "NORMAL_UPDATED";  break;
9812       case UPD_NORMAL_CONV: dump << "UPD_NORMAL_CONV"; break;
9813       case MARKED:          dump << "MARKED";          break;
9814       case MULTI_NORMAL:    dump << "MULTI_NORMAL";    break;
9815       case NEAR_BOUNDARY:   dump << "NEAR_BOUNDARY";   break;
9816       case SMOOTHED_C1:     dump << "SMOOTHED_C1";     break;
9817       case DISTORTED:       dump << "DISTORTED";       break;
9818       case RISKY_SWOL:      dump << "RISKY_SWOL";      break;
9819       case SHRUNK:          dump << "SHRUNK";          break;
9820       case UNUSED_FLAG:     dump << "UNUSED_FLAG";     break;
9821       }
9822       dump << " ";
9823     }
9824   cout << dump << endl;
9825   return dump;
9826 }
9827
9828
9829 //================================================================================
9830 /*!
9831  * \brief Create layers of prisms
9832  */
9833 //================================================================================
9834
9835 bool _ViscousBuilder::refine(_SolidData& data)
9836 {
9837   SMESH_MesherHelper& helper = data.GetHelper();
9838   helper.SetElementsOnShape(false);
9839
9840   Handle(Geom_Curve) curve;
9841   Handle(ShapeAnalysis_Surface) surface;
9842   TopoDS_Edge geomEdge;
9843   TopoDS_Face geomFace;
9844   TopLoc_Location loc;
9845   double f,l, u = 0;
9846   gp_XY uv;
9847   vector< gp_XYZ > pos3D;
9848   bool isOnEdge, isTooConvexFace = false;
9849   TGeomID prevBaseId = -1;
9850   TNode2Edge* n2eMap = 0;
9851   TNode2Edge::iterator n2e;
9852
9853   // Create intermediate nodes on each _LayerEdge
9854
9855   for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
9856   {
9857     _EdgesOnShape& eos = data._edgesOnShape[iS];
9858     if ( eos._edges.empty() ) continue;
9859
9860     if ( eos._edges[0]->_nodes.size() < 2 )
9861       continue; // on _noShrinkShapes
9862
9863     // get data of a shrink shape
9864     isOnEdge = false;
9865     geomEdge.Nullify(); geomFace.Nullify();
9866     curve.Nullify(); surface.Nullify();
9867     if ( !eos._sWOL.IsNull() )
9868     {
9869       isOnEdge = ( eos.SWOLType() == TopAbs_EDGE );
9870       if ( isOnEdge )
9871       {
9872         geomEdge = TopoDS::Edge( eos._sWOL );
9873         curve    = BRep_Tool::Curve( geomEdge, loc, f,l);
9874       }
9875       else
9876       {
9877         geomFace = TopoDS::Face( eos._sWOL );
9878         surface  = helper.GetSurface( geomFace );
9879       }
9880     }
9881     else if ( eos.ShapeType() == TopAbs_FACE && eos._toSmooth )
9882     {
9883       geomFace = TopoDS::Face( eos._shape );
9884       surface  = helper.GetSurface( geomFace );
9885       // propagate _toSmooth back to _eosC1, which was unset in findShapesToSmooth()
9886       for ( size_t i = 0; i < eos._eosC1.size(); ++i )
9887         eos._eosC1[ i ]->_toSmooth = true;
9888
9889       isTooConvexFace = false;
9890       if ( _ConvexFace* cf = data.GetConvexFace( eos._shapeID ))
9891         isTooConvexFace = cf->_isTooCurved;
9892     }
9893
9894     vector< double > segLen;
9895     for ( size_t i = 0; i < eos._edges.size(); ++i )
9896     {
9897       _LayerEdge& edge = *eos._edges[i];
9898       if ( edge._pos.size() < 2 )
9899         continue;
9900
9901       // get accumulated length of segments
9902       segLen.resize( edge._pos.size() );
9903       segLen[0] = 0.0;
9904       if ( eos._sWOL.IsNull() )
9905       {
9906         bool useNormal = true;
9907         bool    usePos = false;
9908         bool  smoothed = false;
9909         double   preci = 0.1 * edge._len;
9910         if ( eos._toSmooth && edge._pos.size() > 2 )
9911         {
9912           smoothed = edge.GetSmoothedPos( preci );
9913         }
9914         if ( smoothed )
9915         {
9916           if ( !surface.IsNull() && !isTooConvexFace ) // edge smoothed on FACE
9917           {
9918             useNormal = usePos = false;
9919             gp_Pnt2d uv = helper.GetNodeUV( geomFace, edge._nodes[0] );
9920             for ( size_t j = 1; j < edge._pos.size() && !useNormal; ++j )
9921             {
9922               uv = surface->NextValueOfUV( uv, edge._pos[j], preci );
9923               if ( surface->Gap() < 2. * edge._len )
9924                 segLen[j] = surface->Gap();
9925               else
9926                 useNormal = true;
9927             }
9928           }
9929         }
9930         else if ( !edge.Is( _LayerEdge::NORMAL_UPDATED ))
9931         {
9932 #ifndef __NODES_AT_POS
9933           useNormal = usePos = false;
9934           edge._pos[1] = edge._pos.back();
9935           edge._pos.resize( 2 );
9936           segLen.resize( 2 );
9937           segLen[ 1 ] = edge._len;
9938 #endif
9939         }
9940         if ( useNormal && edge.Is( _LayerEdge::NORMAL_UPDATED ))
9941         {
9942           useNormal = usePos = false;
9943           _LayerEdge tmpEdge; // get original _normal
9944           tmpEdge._nodes.push_back( edge._nodes[0] );
9945           if ( !setEdgeData( tmpEdge, eos, helper, data ))
9946             usePos = true;
9947           else
9948             for ( size_t j = 1; j < edge._pos.size(); ++j )
9949               segLen[j] = ( edge._pos[j] - edge._pos[0] ) * tmpEdge._normal;
9950         }
9951         if ( useNormal )
9952         {
9953           for ( size_t j = 1; j < edge._pos.size(); ++j )
9954             segLen[j] = ( edge._pos[j] - edge._pos[0] ) * edge._normal;
9955         }
9956         if ( usePos )
9957         {
9958           for ( size_t j = 1; j < edge._pos.size(); ++j )
9959             segLen[j] = segLen[j-1] + ( edge._pos[j-1] - edge._pos[j] ).Modulus();
9960         }
9961         else
9962         {
9963           bool swapped = ( edge._pos.size() > 2 );
9964           while ( swapped )
9965           {
9966             swapped = false;
9967             for ( size_t j = 1; j < edge._pos.size()-1; ++j )
9968               if ( segLen[j] > segLen.back() )
9969               {
9970                 segLen.erase( segLen.begin() + j );
9971                 edge._pos.erase( edge._pos.begin() + j );
9972                 --j;
9973               }
9974               else if ( segLen[j] < segLen[j-1] )
9975               {
9976                 std::swap( segLen[j], segLen[j-1] );
9977                 std::swap( edge._pos[j], edge._pos[j-1] );
9978                 swapped = true;
9979               }
9980           }
9981         }
9982         // smooth a path formed by edge._pos
9983 #ifndef __NODES_AT_POS
9984         if (( smoothed ) /*&&
9985             ( eos.ShapeType() == TopAbs_FACE || edge.Is( _LayerEdge::SMOOTHED_C1 ))*/)
9986           edge.SmoothPos( segLen, preci );
9987 #endif
9988       }
9989       else if ( eos._isRegularSWOL ) // usual SWOL
9990       {
9991         if ( edge.Is( _LayerEdge::SMOOTHED ))
9992         {
9993           SMESH_NodeXYZ p0( edge._nodes[0] );
9994           for ( size_t j = 1; j < edge._pos.size(); ++j )
9995           {
9996             gp_XYZ pj = surface->Value( edge._pos[j].X(), edge._pos[j].Y() ).XYZ();
9997             segLen[j] = ( pj - p0 ) * edge._normal;
9998           }
9999         }
10000         else
10001         {
10002           for ( size_t j = 1; j < edge._pos.size(); ++j )
10003             segLen[j] = segLen[j-1] + (edge._pos[j-1] - edge._pos[j] ).Modulus();
10004         }
10005       }
10006       else if ( !surface.IsNull() ) // SWOL surface with singularities
10007       {
10008         pos3D.resize( edge._pos.size() );
10009         for ( size_t j = 0; j < edge._pos.size(); ++j )
10010           pos3D[j] = surface->Value( edge._pos[j].X(), edge._pos[j].Y() ).XYZ();
10011
10012         for ( size_t j = 1; j < edge._pos.size(); ++j )
10013           segLen[j] = segLen[j-1] + ( pos3D[j-1] - pos3D[j] ).Modulus();
10014       }
10015
10016       // allocate memory for new nodes if it is not yet refined
10017       const SMDS_MeshNode* tgtNode = edge._nodes.back();
10018       if ( edge._nodes.size() == 2 )
10019       {
10020 #ifdef __NODES_AT_POS
10021         int nbNodes = edge._pos.size();
10022 #else
10023         int nbNodes = eos._hyp.GetNumberLayers() + 1;
10024 #endif
10025         edge._nodes.resize( nbNodes, 0 );
10026         edge._nodes[1] = 0;
10027         edge._nodes.back() = tgtNode;
10028       }
10029       // restore shapePos of the last node by already treated _LayerEdge of another _SolidData
10030       const TGeomID baseShapeId = edge._nodes[0]->getshapeId();
10031       if ( baseShapeId != prevBaseId )
10032       {
10033         map< TGeomID, TNode2Edge* >::iterator s2ne = data._s2neMap.find( baseShapeId );
10034         n2eMap = ( s2ne == data._s2neMap.end() ) ? 0 : s2ne->second;
10035         prevBaseId = baseShapeId;
10036       }
10037       _LayerEdge* edgeOnSameNode = 0;
10038       bool        useExistingPos = false;
10039       if ( n2eMap && (( n2e = n2eMap->find( edge._nodes[0] )) != n2eMap->end() ))
10040       {
10041         edgeOnSameNode = n2e->second;
10042         useExistingPos = ( edgeOnSameNode->_len < edge._len );
10043         const gp_XYZ& otherTgtPos = edgeOnSameNode->_pos.back();
10044         SMDS_PositionPtr  lastPos = tgtNode->GetPosition();
10045         if ( isOnEdge )
10046         {
10047           SMDS_EdgePositionPtr epos = lastPos;
10048           epos->SetUParameter( otherTgtPos.X() );
10049         }
10050         else
10051         {
10052           SMDS_FacePositionPtr fpos = lastPos;
10053           fpos->SetUParameter( otherTgtPos.X() );
10054           fpos->SetVParameter( otherTgtPos.Y() );
10055         }
10056       }
10057       // calculate height of the first layer
10058       double h0;
10059       const double T = segLen.back(); //data._hyp.GetTotalThickness();
10060       const double f = eos._hyp.GetStretchFactor();
10061       const int    N = eos._hyp.GetNumberLayers();
10062       const double fPowN = pow( f, N );
10063       if ( fPowN - 1 <= numeric_limits<double>::min() )
10064         h0 = T / N;
10065       else
10066         h0 = T * ( f - 1 )/( fPowN - 1 );
10067
10068       const double zeroLen = std::numeric_limits<double>::min();
10069
10070       // create intermediate nodes
10071       double hSum = 0, hi = h0/f;
10072       size_t iSeg = 1;
10073       for ( size_t iStep = 1; iStep < edge._nodes.size(); ++iStep )
10074       {
10075         // compute an intermediate position
10076         hi *= f;
10077         hSum += hi;
10078         while ( hSum > segLen[iSeg] && iSeg < segLen.size()-1 )
10079           ++iSeg;
10080         int iPrevSeg = iSeg-1;
10081         while ( fabs( segLen[iPrevSeg] - segLen[iSeg]) <= zeroLen && iPrevSeg > 0 )
10082           --iPrevSeg;
10083         double   r = ( segLen[iSeg] - hSum ) / ( segLen[iSeg] - segLen[iPrevSeg] );
10084         gp_Pnt pos = r * edge._pos[iPrevSeg] + (1-r) * edge._pos[iSeg];
10085 #ifdef __NODES_AT_POS
10086         pos = edge._pos[ iStep ];
10087 #endif
10088         SMDS_MeshNode*& node = const_cast< SMDS_MeshNode*& >( edge._nodes[ iStep ]);
10089         if ( !eos._sWOL.IsNull() )
10090         {
10091           // compute XYZ by parameters <pos>
10092           if ( isOnEdge )
10093           {
10094             u = pos.X();
10095             if ( !node )
10096               pos = curve->Value( u ).Transformed(loc);
10097           }
10098           else if ( eos._isRegularSWOL )
10099           {
10100             uv.SetCoord( pos.X(), pos.Y() );
10101             if ( !node )
10102               pos = surface->Value( pos.X(), pos.Y() );
10103           }
10104           else
10105           {
10106             uv.SetCoord( pos.X(), pos.Y() );
10107             gp_Pnt p = r * pos3D[ iPrevSeg ] + (1-r) * pos3D[ iSeg ];
10108             uv = surface->NextValueOfUV( uv, p, BRep_Tool::Tolerance( geomFace )).XY();
10109             if ( !node )
10110               pos = surface->Value( uv );
10111           }
10112         }
10113         // create or update the node
10114         if ( !node )
10115         {
10116           node = helper.AddNode( pos.X(), pos.Y(), pos.Z());
10117           if ( !eos._sWOL.IsNull() )
10118           {
10119             if ( isOnEdge )
10120               getMeshDS()->SetNodeOnEdge( node, geomEdge, u );
10121             else
10122               getMeshDS()->SetNodeOnFace( node, geomFace, uv.X(), uv.Y() );
10123           }
10124           else
10125           {
10126             getMeshDS()->SetNodeInVolume( node, helper.GetSubShapeID() );
10127           }
10128         }
10129         else
10130         {
10131           if ( !eos._sWOL.IsNull() )
10132           {
10133             // make average pos from new and current parameters
10134             if ( isOnEdge )
10135             {
10136               //u = 0.5 * ( u + helper.GetNodeU( geomEdge, node ));
10137               if ( useExistingPos )
10138                 u = helper.GetNodeU( geomEdge, node );
10139               pos = curve->Value( u ).Transformed(loc);
10140
10141               SMDS_EdgePositionPtr epos = node->GetPosition();
10142               epos->SetUParameter( u );
10143             }
10144             else
10145             {
10146               //uv = 0.5 * ( uv + helper.GetNodeUV( geomFace, node ));
10147               if ( useExistingPos )
10148                 uv = helper.GetNodeUV( geomFace, node );
10149               pos = surface->Value( uv );
10150
10151               SMDS_FacePositionPtr fpos = node->GetPosition();
10152               fpos->SetUParameter( uv.X() );
10153               fpos->SetVParameter( uv.Y() );
10154             }
10155           }
10156           node->setXYZ( pos.X(), pos.Y(), pos.Z() );
10157         }
10158       } // loop on edge._nodes
10159
10160       if ( !eos._sWOL.IsNull() ) // prepare for shrink()
10161       {
10162         if ( isOnEdge )
10163           edge._pos.back().SetCoord( u, 0,0);
10164         else
10165           edge._pos.back().SetCoord( uv.X(), uv.Y() ,0);
10166
10167         if ( edgeOnSameNode )
10168           edgeOnSameNode->_pos.back() = edge._pos.back();
10169       }
10170
10171     } // loop on eos._edges to create nodes
10172
10173
10174     if ( !getMeshDS()->IsEmbeddedMode() )
10175       // Log node movement
10176       for ( size_t i = 0; i < eos._edges.size(); ++i )
10177       {
10178         SMESH_TNodeXYZ p ( eos._edges[i]->_nodes.back() );
10179         getMeshDS()->MoveNode( p._node, p.X(), p.Y(), p.Z() );
10180       }
10181   }
10182
10183
10184   // Create volumes
10185
10186   helper.SetElementsOnShape(true);
10187
10188   vector< vector<const SMDS_MeshNode*>* > nnVec;
10189   set< vector<const SMDS_MeshNode*>* >    nnSet;
10190   set< int >                       degenEdgeInd;
10191   vector<const SMDS_MeshElement*>     degenVols;
10192
10193   TopExp_Explorer exp( data._solid, TopAbs_FACE );
10194   for ( ; exp.More(); exp.Next() )
10195   {
10196     const TGeomID faceID = getMeshDS()->ShapeToIndex( exp.Current() );
10197     if ( data._ignoreFaceIds.count( faceID ))
10198       continue;
10199     _EdgesOnShape*    eos = data.GetShapeEdges( faceID );
10200     SMDS_MeshGroup* group = StdMeshers_ViscousLayers::CreateGroup( eos->_hyp.GetGroupName(),
10201                                                                    *helper.GetMesh(),
10202                                                                    SMDSAbs_Volume );
10203     std::vector< const SMDS_MeshElement* > vols;
10204     const bool isReversedFace = data._reversedFaceIds.count( faceID );
10205     SMESHDS_SubMesh*    fSubM = getMeshDS()->MeshElements( exp.Current() );
10206     SMDS_ElemIteratorPtr  fIt = fSubM->GetElements();
10207     while ( fIt->more() )
10208     {
10209       const SMDS_MeshElement* face = fIt->next();
10210       const int            nbNodes = face->NbCornerNodes();
10211       nnVec.resize( nbNodes );
10212       nnSet.clear();
10213       degenEdgeInd.clear();
10214       size_t maxZ = 0, minZ = std::numeric_limits<size_t>::max();
10215       SMDS_NodeIteratorPtr nIt = face->nodeIterator();
10216       for ( int iN = 0; iN < nbNodes; ++iN )
10217       {
10218         const SMDS_MeshNode* n = nIt->next();
10219         _LayerEdge*       edge = data._n2eMap[ n ];
10220         const int i = isReversedFace ? nbNodes-1-iN : iN;
10221         nnVec[ i ] = & edge->_nodes;
10222         maxZ = std::max( maxZ, nnVec[ i ]->size() );
10223         minZ = std::min( minZ, nnVec[ i ]->size() );
10224
10225         if ( helper.HasDegeneratedEdges() )
10226           nnSet.insert( nnVec[ i ]);
10227       }
10228
10229       if ( maxZ == 0 )
10230         continue;
10231       if ( 0 < nnSet.size() && nnSet.size() < 3 )
10232         continue;
10233
10234       vols.clear();
10235       const SMDS_MeshElement* vol;
10236
10237       switch ( nbNodes )
10238       {
10239       case 3: // TRIA
10240       {
10241         // PENTA
10242         for ( size_t iZ = 1; iZ < minZ; ++iZ )
10243         {
10244           vol = helper.AddVolume( (*nnVec[0])[iZ-1], (*nnVec[1])[iZ-1], (*nnVec[2])[iZ-1],
10245                                   (*nnVec[0])[iZ],   (*nnVec[1])[iZ],   (*nnVec[2])[iZ]);
10246           vols.push_back( vol );
10247         }
10248
10249         for ( size_t iZ = minZ; iZ < maxZ; ++iZ )
10250         {
10251           for ( int iN = 0; iN < nbNodes; ++iN )
10252             if ( nnVec[ iN ]->size() < iZ+1 )
10253               degenEdgeInd.insert( iN );
10254
10255           if ( degenEdgeInd.size() == 1 )  // PYRAM
10256           {
10257             int i2 = *degenEdgeInd.begin();
10258             int i0 = helper.WrapIndex( i2 - 1, nbNodes );
10259             int i1 = helper.WrapIndex( i2 + 1, nbNodes );
10260             vol = helper.AddVolume( (*nnVec[i0])[iZ-1], (*nnVec[i1])[iZ-1],
10261                                     (*nnVec[i1])[iZ  ], (*nnVec[i0])[iZ  ], (*nnVec[i2]).back());
10262             vols.push_back( vol );
10263           }
10264           else  // TETRA
10265           {
10266             int i3 = !degenEdgeInd.count(0) ? 0 : !degenEdgeInd.count(1) ? 1 : 2;
10267             vol = helper.AddVolume( (*nnVec[  0 ])[ i3 == 0 ? iZ-1 : nnVec[0]->size()-1 ],
10268                                     (*nnVec[  1 ])[ i3 == 1 ? iZ-1 : nnVec[1]->size()-1 ],
10269                                     (*nnVec[  2 ])[ i3 == 2 ? iZ-1 : nnVec[2]->size()-1 ],
10270                                     (*nnVec[ i3 ])[ iZ ]);
10271             vols.push_back( vol );
10272           }
10273         }
10274         break; // TRIA
10275       }
10276       case 4: // QUAD
10277       {
10278         // HEX
10279         for ( size_t iZ = 1; iZ < minZ; ++iZ )
10280         {
10281           vol = helper.AddVolume( (*nnVec[0])[iZ-1], (*nnVec[1])[iZ-1],
10282                                   (*nnVec[2])[iZ-1], (*nnVec[3])[iZ-1],
10283                                   (*nnVec[0])[iZ],   (*nnVec[1])[iZ],
10284                                   (*nnVec[2])[iZ],   (*nnVec[3])[iZ]);
10285           vols.push_back( vol );
10286         }
10287
10288         for ( size_t iZ = minZ; iZ < maxZ; ++iZ )
10289         {
10290           for ( int iN = 0; iN < nbNodes; ++iN )
10291             if ( nnVec[ iN ]->size() < iZ+1 )
10292               degenEdgeInd.insert( iN );
10293
10294           switch ( degenEdgeInd.size() )
10295           {
10296           case 2: // PENTA
10297           {
10298             int i2 = *degenEdgeInd.begin();
10299             int i3 = *degenEdgeInd.rbegin();
10300             bool ok = ( i3 - i2 == 1 );
10301             if ( i2 == 0 && i3 == 3 ) { i2 = 3; i3 = 0; ok = true; }
10302             int i0 = helper.WrapIndex( i3 + 1, nbNodes );
10303             int i1 = helper.WrapIndex( i0 + 1, nbNodes );
10304
10305             vol = helper.AddVolume( nnVec[i3]->back(), (*nnVec[i0])[iZ], (*nnVec[i0])[iZ-1],
10306                                     nnVec[i2]->back(), (*nnVec[i1])[iZ], (*nnVec[i1])[iZ-1]);
10307             vols.push_back( vol );
10308             if ( !ok && vol )
10309               degenVols.push_back( vol );
10310           }
10311           break;
10312
10313           default: // degen HEX
10314           {
10315             vol = helper.AddVolume( nnVec[0]->size() > iZ-1 ? (*nnVec[0])[iZ-1] : nnVec[0]->back(),
10316                                     nnVec[1]->size() > iZ-1 ? (*nnVec[1])[iZ-1] : nnVec[1]->back(),
10317                                     nnVec[2]->size() > iZ-1 ? (*nnVec[2])[iZ-1] : nnVec[2]->back(),
10318                                     nnVec[3]->size() > iZ-1 ? (*nnVec[3])[iZ-1] : nnVec[3]->back(),
10319                                     nnVec[0]->size() > iZ   ? (*nnVec[0])[iZ]   : nnVec[0]->back(),
10320                                     nnVec[1]->size() > iZ   ? (*nnVec[1])[iZ]   : nnVec[1]->back(),
10321                                     nnVec[2]->size() > iZ   ? (*nnVec[2])[iZ]   : nnVec[2]->back(),
10322                                     nnVec[3]->size() > iZ   ? (*nnVec[3])[iZ]   : nnVec[3]->back());
10323             vols.push_back( vol );
10324             degenVols.push_back( vol );
10325           }
10326           }
10327         }
10328         break; // HEX
10329       }
10330       default:
10331         return error("Not supported type of element", data._index);
10332
10333       } // switch ( nbNodes )
10334
10335       if ( group )
10336         for ( size_t i = 0; i < vols.size(); ++i )
10337           group->Add( vols[ i ]);
10338
10339     } // while ( fIt->more() )
10340   } // loop on FACEs
10341
10342   if ( !degenVols.empty() )
10343   {
10344     SMESH_ComputeErrorPtr& err = _mesh->GetSubMesh( data._solid )->GetComputeError();
10345     if ( !err || err->IsOK() )
10346     {
10347       SMESH_BadInputElements* badElems =
10348         new SMESH_BadInputElements( getMeshDS(), COMPERR_WARNING, "Bad quality volumes created" );
10349       badElems->myBadElements.insert( badElems->myBadElements.end(),
10350                                       degenVols.begin(),degenVols.end() );
10351       err.reset( badElems );
10352     }
10353   }
10354
10355   return true;
10356 }
10357
10358 //================================================================================
10359 /*!
10360  * \brief Shrink 2D mesh on faces to let space for inflated layers
10361  */
10362 //================================================================================
10363
10364 bool _ViscousBuilder::shrink(_SolidData& theData)
10365 {
10366   // make map of (ids of FACEs to shrink mesh on) to (list of _SolidData containing
10367   // _LayerEdge's inflated along FACE or EDGE)
10368   map< TGeomID, list< _SolidData* > > f2sdMap;
10369   for ( size_t i = 0 ; i < _sdVec.size(); ++i )
10370   {
10371     _SolidData& data = _sdVec[i];
10372     map< TGeomID, TopoDS_Shape >::iterator s2s = data._shrinkShape2Shape.begin();
10373     for (; s2s != data._shrinkShape2Shape.end(); ++s2s )
10374       if ( s2s->second.ShapeType() == TopAbs_FACE && !_shrinkedFaces.Contains( s2s->second ))
10375       {
10376         f2sdMap[ getMeshDS()->ShapeToIndex( s2s->second )].push_back( &data );
10377
10378         // Put mesh faces on the shrinked FACE to the proxy sub-mesh to avoid
10379         // usage of mesh faces made in addBoundaryElements() by the 3D algo or
10380         // by StdMeshers_QuadToTriaAdaptor
10381         if ( SMESHDS_SubMesh* smDS = getMeshDS()->MeshElements( s2s->second ))
10382         {
10383           SMESH_ProxyMesh::SubMesh* proxySub =
10384             data._proxyMesh->getFaceSubM( TopoDS::Face( s2s->second ), /*create=*/true);
10385           if ( proxySub->NbElements() == 0 )
10386           {
10387             SMDS_ElemIteratorPtr fIt = smDS->GetElements();
10388             while ( fIt->more() )
10389             {
10390               const SMDS_MeshElement* f = fIt->next();
10391               // as a result 3D algo will use elements from proxySub and not from smDS
10392               proxySub->AddElement( f );
10393               f->setIsMarked( true );
10394
10395               // Mark nodes on the FACE to discriminate them from nodes
10396               // added by addBoundaryElements(); marked nodes are to be smoothed while shrink()
10397               for ( int iN = 0, nbN = f->NbNodes(); iN < nbN; ++iN )
10398               {
10399                 const SMDS_MeshNode* n = f->GetNode( iN );
10400                 if ( n->GetPosition()->GetDim() == 2 )
10401                   n->setIsMarked( true );
10402               }
10403             }
10404           }
10405         }
10406       }
10407   }
10408
10409   SMESH_MesherHelper helper( *_mesh );
10410   helper.ToFixNodeParameters( true );
10411
10412   // EDGEs to shrink
10413   map< TGeomID, _Shrinker1D > e2shrMap;
10414   vector< _EdgesOnShape* > subEOS;
10415   vector< _LayerEdge* > lEdges;
10416
10417   // loop on FACEs to shrink mesh on
10418   map< TGeomID, list< _SolidData* > >::iterator f2sd = f2sdMap.begin();
10419   for ( ; f2sd != f2sdMap.end(); ++f2sd )
10420   {
10421     list< _SolidData* > & dataList = f2sd->second;
10422     if ( dataList.front()->_n2eMap.empty() ||
10423          dataList.back() ->_n2eMap.empty() )
10424       continue; // not yet computed
10425     if ( dataList.front() != &theData &&
10426          dataList.back()  != &theData )
10427       continue;
10428
10429     _SolidData&      data = *dataList.front();
10430     _SolidData*     data2 = dataList.size() > 1 ? dataList.back() : 0;
10431     const TopoDS_Face&  F = TopoDS::Face( getMeshDS()->IndexToShape( f2sd->first ));
10432     SMESH_subMesh*     sm = _mesh->GetSubMesh( F );
10433     SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
10434
10435     Handle(Geom_Surface) surface = BRep_Tool::Surface( F );
10436
10437     _shrinkedFaces.Add( F );
10438     helper.SetSubShape( F );
10439
10440     // ===========================
10441     // Prepare data for shrinking
10442     // ===========================
10443
10444     // Collect nodes to smooth (they are marked at the beginning of this method)
10445     vector < const SMDS_MeshNode* > smoothNodes;
10446     {
10447       SMDS_NodeIteratorPtr nIt = smDS->GetNodes();
10448       while ( nIt->more() )
10449       {
10450         const SMDS_MeshNode* n = nIt->next();
10451         if ( n->isMarked() )
10452           smoothNodes.push_back( n );
10453       }
10454     }
10455     // Find out face orientation
10456     double refSign = 1;
10457     const set<TGeomID> ignoreShapes;
10458     bool isOkUV;
10459     if ( !smoothNodes.empty() )
10460     {
10461       vector<_Simplex> simplices;
10462       _Simplex::GetSimplices( smoothNodes[0], simplices, ignoreShapes );
10463       helper.GetNodeUV( F, simplices[0]._nPrev, 0, &isOkUV ); // fix UV of simplex nodes
10464       helper.GetNodeUV( F, simplices[0]._nNext, 0, &isOkUV );
10465       gp_XY uv = helper.GetNodeUV( F, smoothNodes[0], 0, &isOkUV );
10466       if ( !simplices[0].IsForward(uv, smoothNodes[0], F, helper, refSign ))
10467         refSign = -1;
10468     }
10469
10470     // Find _LayerEdge's inflated along F
10471     subEOS.clear();
10472     lEdges.clear();
10473     {
10474       SMESH_subMeshIteratorPtr subIt = sm->getDependsOnIterator(/*includeSelf=*/false,
10475                                                                 /*complexFirst=*/true); //!!!
10476       while ( subIt->more() )
10477       {
10478         const TGeomID subID = subIt->next()->GetId();
10479         if ( data._noShrinkShapes.count( subID ))
10480           continue;
10481         _EdgesOnShape* eos = data.GetShapeEdges( subID );
10482         if ( !eos || eos->_sWOL.IsNull() )
10483           if ( data2 ) // check in adjacent SOLID
10484           {
10485             eos = data2->GetShapeEdges( subID );
10486             if ( !eos || eos->_sWOL.IsNull() )
10487               continue;
10488           }
10489         subEOS.push_back( eos );
10490
10491         for ( size_t i = 0; i < eos->_edges.size(); ++i )
10492         {
10493           lEdges.push_back( eos->_edges[ i ] );
10494           prepareEdgeToShrink( *eos->_edges[ i ], *eos, helper, smDS );
10495         }
10496       }
10497     }
10498
10499     dumpFunction(SMESH_Comment("beforeShrinkFace")<<f2sd->first); // debug
10500     SMDS_ElemIteratorPtr fIt = smDS->GetElements();
10501     while ( fIt->more() )
10502       if ( const SMDS_MeshElement* f = fIt->next() )
10503         dumpChangeNodes( f );
10504     dumpFunctionEnd();
10505
10506     // Replace source nodes by target nodes in mesh faces to shrink
10507     dumpFunction(SMESH_Comment("replNodesOnFace")<<f2sd->first); // debug
10508     const SMDS_MeshNode* nodes[20];
10509     for ( size_t iS = 0; iS < subEOS.size(); ++iS )
10510     {
10511       _EdgesOnShape& eos = * subEOS[ iS ];
10512       for ( size_t i = 0; i < eos._edges.size(); ++i )
10513       {
10514         _LayerEdge& edge = *eos._edges[i];
10515         const SMDS_MeshNode* srcNode = edge._nodes[0];
10516         const SMDS_MeshNode* tgtNode = edge._nodes.back();
10517         SMDS_ElemIteratorPtr fIt = srcNode->GetInverseElementIterator(SMDSAbs_Face);
10518         while ( fIt->more() )
10519         {
10520           const SMDS_MeshElement* f = fIt->next();
10521           if ( !smDS->Contains( f ) || !f->isMarked() )
10522             continue;
10523           SMDS_NodeIteratorPtr nIt = f->nodeIterator();
10524           for ( int iN = 0; nIt->more(); ++iN )
10525           {
10526             const SMDS_MeshNode* n = nIt->next();
10527             nodes[iN] = ( n == srcNode ? tgtNode : n );
10528           }
10529           helper.GetMeshDS()->ChangeElementNodes( f, nodes, f->NbNodes() );
10530           dumpChangeNodes( f );
10531         }
10532       }
10533     }
10534     dumpFunctionEnd();
10535
10536     // find out if a FACE is concave
10537     const bool isConcaveFace = isConcave( F, helper );
10538
10539     // Create _SmoothNode's on face F
10540     vector< _SmoothNode > nodesToSmooth( smoothNodes.size() );
10541     {
10542       dumpFunction(SMESH_Comment("fixUVOnFace")<<f2sd->first); // debug
10543       const bool sortSimplices = isConcaveFace;
10544       for ( size_t i = 0; i < smoothNodes.size(); ++i )
10545       {
10546         const SMDS_MeshNode* n = smoothNodes[i];
10547         nodesToSmooth[ i ]._node = n;
10548         // src nodes must be already replaced by tgt nodes to have tgt nodes in _simplices
10549         _Simplex::GetSimplices( n, nodesToSmooth[ i ]._simplices, ignoreShapes, 0, sortSimplices);
10550         // fix up incorrect uv of nodes on the FACE
10551         helper.GetNodeUV( F, n, 0, &isOkUV);
10552         dumpMove( n );
10553       }
10554       dumpFunctionEnd();
10555     }
10556     //if ( nodesToSmooth.empty() ) continue;
10557
10558     // Find EDGE's to shrink and set simpices to LayerEdge's
10559     set< _Shrinker1D* > eShri1D;
10560     {
10561       for ( size_t iS = 0; iS < subEOS.size(); ++iS )
10562       {
10563         _EdgesOnShape& eos = * subEOS[ iS ];
10564         if ( eos.SWOLType() == TopAbs_EDGE )
10565         {
10566           SMESH_subMesh* edgeSM = _mesh->GetSubMesh( eos._sWOL );
10567           _Shrinker1D& shrinker  = e2shrMap[ edgeSM->GetId() ];
10568           eShri1D.insert( & shrinker );
10569           shrinker.AddEdge( eos._edges[0], eos, helper );
10570           VISCOUS_3D::ToClearSubWithMain( edgeSM, data._solid );
10571           // restore params of nodes on EDGE if the EDGE has been already
10572           // shrinked while shrinking other FACE
10573           shrinker.RestoreParams();
10574         }
10575         for ( size_t i = 0; i < eos._edges.size(); ++i )
10576         {
10577           _LayerEdge& edge = * eos._edges[i];
10578           _Simplex::GetSimplices( /*tgtNode=*/edge._nodes.back(), edge._simplices, ignoreShapes );
10579
10580           // additionally mark tgt node; only marked nodes will be used in SetNewLength2d()
10581           // not-marked nodes are those added by refine()
10582           edge._nodes.back()->setIsMarked( true );
10583         }
10584       }
10585     }
10586
10587     bool toFixTria = false; // to improve quality of trias by diagonal swap
10588     if ( isConcaveFace )
10589     {
10590       const bool hasTria = _mesh->NbTriangles(), hasQuad = _mesh->NbQuadrangles();
10591       if ( hasTria != hasQuad ) {
10592         toFixTria = hasTria;
10593       }
10594       else {
10595         set<int> nbNodesSet;
10596         SMDS_ElemIteratorPtr fIt = smDS->GetElements();
10597         while ( fIt->more() && nbNodesSet.size() < 2 )
10598           nbNodesSet.insert( fIt->next()->NbCornerNodes() );
10599         toFixTria = ( *nbNodesSet.begin() == 3 );
10600       }
10601     }
10602
10603     // ==================
10604     // Perform shrinking
10605     // ==================
10606
10607     bool shrinked = true;
10608     int nbBad, shriStep=0, smooStep=0;
10609     _SmoothNode::SmoothType smoothType
10610       = isConcaveFace ? _SmoothNode::ANGULAR : _SmoothNode::LAPLACIAN;
10611     SMESH_Comment errMsg;
10612     while ( shrinked )
10613     {
10614       shriStep++;
10615       // Move boundary nodes (actually just set new UV)
10616       // -----------------------------------------------
10617       dumpFunction(SMESH_Comment("moveBoundaryOnF")<<f2sd->first<<"_st"<<shriStep ); // debug
10618       shrinked = false;
10619       for ( size_t iS = 0; iS < subEOS.size(); ++iS )
10620       {
10621         _EdgesOnShape& eos = * subEOS[ iS ];
10622         for ( size_t i = 0; i < eos._edges.size(); ++i )
10623         {
10624           shrinked |= eos._edges[i]->SetNewLength2d( surface, F, eos, helper );
10625         }
10626       }
10627       dumpFunctionEnd();
10628
10629       // Move nodes on EDGE's
10630       // (XYZ is set as soon as a needed length reached in SetNewLength2d())
10631       set< _Shrinker1D* >::iterator shr = eShri1D.begin();
10632       for ( ; shr != eShri1D.end(); ++shr )
10633         (*shr)->Compute( /*set3D=*/false, helper );
10634
10635       // Smoothing in 2D
10636       // -----------------
10637       int nbNoImpSteps = 0;
10638       bool       moved = true;
10639       nbBad = 1;
10640       while (( nbNoImpSteps < 5 && nbBad > 0) && moved)
10641       {
10642         dumpFunction(SMESH_Comment("shrinkFace")<<f2sd->first<<"_st"<<++smooStep); // debug
10643
10644         int oldBadNb = nbBad;
10645         nbBad = 0;
10646         moved = false;
10647         // '% 5' minimizes NB FUNCTIONS on viscous_layers_00/B2 case
10648         _SmoothNode::SmoothType smooTy = ( smooStep % 5 ) ? smoothType : _SmoothNode::LAPLACIAN;
10649         for ( size_t i = 0; i < nodesToSmooth.size(); ++i )
10650         {
10651           moved |= nodesToSmooth[i].Smooth( nbBad, surface, helper, refSign,
10652                                             smooTy, /*set3D=*/isConcaveFace);
10653         }
10654         if ( nbBad < oldBadNb )
10655           nbNoImpSteps = 0;
10656         else
10657           nbNoImpSteps++;
10658
10659         dumpFunctionEnd();
10660       }
10661
10662       errMsg.clear();
10663       if ( nbBad > 0 )
10664         errMsg << "Can't shrink 2D mesh on face " << f2sd->first;
10665       if ( shriStep > 200 )
10666         errMsg << "Infinite loop at shrinking 2D mesh on face " << f2sd->first;
10667       if ( !errMsg.empty() )
10668         break;
10669
10670       // Fix narrow triangles by swapping diagonals
10671       // ---------------------------------------
10672       if ( toFixTria )
10673       {
10674         set<const SMDS_MeshNode*> usedNodes;
10675         fixBadFaces( F, helper, /*is2D=*/true, shriStep, & usedNodes); // swap diagonals
10676
10677         // update working data
10678         set<const SMDS_MeshNode*>::iterator n;
10679         for ( size_t i = 0; i < nodesToSmooth.size() && !usedNodes.empty(); ++i )
10680         {
10681           n = usedNodes.find( nodesToSmooth[ i ]._node );
10682           if ( n != usedNodes.end())
10683           {
10684             _Simplex::GetSimplices( nodesToSmooth[ i ]._node,
10685                                     nodesToSmooth[ i ]._simplices,
10686                                     ignoreShapes, NULL,
10687                                     /*sortSimplices=*/ smoothType == _SmoothNode::ANGULAR );
10688             usedNodes.erase( n );
10689           }
10690         }
10691         for ( size_t i = 0; i < lEdges.size() && !usedNodes.empty(); ++i )
10692         {
10693           n = usedNodes.find( /*tgtNode=*/ lEdges[i]->_nodes.back() );
10694           if ( n != usedNodes.end())
10695           {
10696             _Simplex::GetSimplices( lEdges[i]->_nodes.back(),
10697                                     lEdges[i]->_simplices,
10698                                     ignoreShapes );
10699             usedNodes.erase( n );
10700           }
10701         }
10702       }
10703       // TODO: check effect of this additional smooth
10704       // additional laplacian smooth to increase allowed shrink step
10705       // for ( int st = 1; st; --st )
10706       // {
10707       //   dumpFunction(SMESH_Comment("shrinkFace")<<f2sd->first<<"_st"<<++smooStep); // debug
10708       //   for ( size_t i = 0; i < nodesToSmooth.size(); ++i )
10709       //   {
10710       //     nodesToSmooth[i].Smooth( nbBad,surface,helper,refSign,
10711       //                              _SmoothNode::LAPLACIAN,/*set3D=*/false);
10712       //   }
10713       // }
10714
10715     } // while ( shrinked )
10716
10717     if ( !errMsg.empty() ) // Try to re-compute the shrink FACE
10718     {
10719       debugMsg( "Re-compute FACE " << f2sd->first << " because " << errMsg );
10720
10721       // remove faces
10722       SMESHDS_SubMesh* psm = data._proxyMesh->getFaceSubM( F );
10723       {
10724         vector< const SMDS_MeshElement* > facesToRm;
10725         if ( psm )
10726         {
10727           facesToRm.reserve( psm->NbElements() );
10728           for ( SMDS_ElemIteratorPtr ite = psm->GetElements(); ite->more(); )
10729             facesToRm.push_back( ite->next() );
10730
10731           for ( size_t i = 0 ; i < _sdVec.size(); ++i )
10732             if (( psm = _sdVec[i]._proxyMesh->getFaceSubM( F )))
10733               psm->Clear();
10734         }
10735         for ( size_t i = 0; i < facesToRm.size(); ++i )
10736           getMeshDS()->RemoveFreeElement( facesToRm[i], smDS, /*fromGroups=*/false );
10737       }
10738       // remove nodes
10739       {
10740         TIDSortedNodeSet nodesToKeep; // nodes of _LayerEdge to keep
10741         for ( size_t iS = 0; iS < subEOS.size(); ++iS ) {
10742           for ( size_t i = 0; i < subEOS[iS]->_edges.size(); ++i )
10743             nodesToKeep.insert( ++( subEOS[iS]->_edges[i]->_nodes.begin() ),
10744                                 subEOS[iS]->_edges[i]->_nodes.end() );
10745         }
10746         SMDS_NodeIteratorPtr itn = smDS->GetNodes();
10747         while ( itn->more() ) {
10748           const SMDS_MeshNode* n = itn->next();
10749           if ( !nodesToKeep.count( n ))
10750             getMeshDS()->RemoveFreeNode( n, smDS, /*fromGroups=*/false );
10751         }
10752       }
10753       // restore position and UV of target nodes
10754       gp_Pnt p;
10755       for ( size_t iS = 0; iS < subEOS.size(); ++iS )
10756         for ( size_t i = 0; i < subEOS[iS]->_edges.size(); ++i )
10757         {
10758           _LayerEdge*       edge = subEOS[iS]->_edges[i];
10759           SMDS_MeshNode* tgtNode = const_cast< SMDS_MeshNode*& >( edge->_nodes.back() );
10760           if ( edge->_pos.empty() ||
10761                edge->Is( _LayerEdge::SHRUNK )) continue;
10762           if ( subEOS[iS]->SWOLType() == TopAbs_FACE )
10763           {
10764             SMDS_FacePositionPtr pos = tgtNode->GetPosition();
10765             pos->SetUParameter( edge->_pos[0].X() );
10766             pos->SetVParameter( edge->_pos[0].Y() );
10767             p = surface->Value( edge->_pos[0].X(), edge->_pos[0].Y() );
10768           }
10769           else
10770           {
10771             SMDS_EdgePositionPtr pos = tgtNode->GetPosition();
10772             pos->SetUParameter( edge->_pos[0].Coord( U_TGT ));
10773             p = BRepAdaptor_Curve( TopoDS::Edge( subEOS[iS]->_sWOL )).Value( pos->GetUParameter() );
10774           }
10775           tgtNode->setXYZ( p.X(), p.Y(), p.Z() );
10776           dumpMove( tgtNode );
10777         }
10778       // shrink EDGE sub-meshes and set proxy sub-meshes
10779       UVPtStructVec uvPtVec;
10780       set< _Shrinker1D* >::iterator shrIt = eShri1D.begin();
10781       for ( shrIt = eShri1D.begin(); shrIt != eShri1D.end(); ++shrIt )
10782       {
10783         _Shrinker1D* shr = (*shrIt);
10784         shr->Compute( /*set3D=*/true, helper );
10785
10786         // set proxy mesh of EDGEs w/o layers
10787         map< double, const SMDS_MeshNode* > nodes;
10788         SMESH_Algo::GetSortedNodesOnEdge( getMeshDS(), shr->GeomEdge(),/*skipMedium=*/true, nodes);
10789         // remove refinement nodes
10790         const SMDS_MeshNode* sn0 = shr->SrcNode(0), *sn1 = shr->SrcNode(1);
10791         const SMDS_MeshNode* tn0 = shr->TgtNode(0), *tn1 = shr->TgtNode(1);
10792         map< double, const SMDS_MeshNode* >::iterator u2n = nodes.begin();
10793         if ( u2n->second == sn0 || u2n->second == sn1 )
10794         {
10795           while ( u2n->second != tn0 && u2n->second != tn1 )
10796             ++u2n;
10797           nodes.erase( nodes.begin(), u2n );
10798         }
10799         u2n = --nodes.end();
10800         if ( u2n->second == sn0 || u2n->second == sn1 )
10801         {
10802           while ( u2n->second != tn0 && u2n->second != tn1 )
10803             --u2n;
10804           nodes.erase( ++u2n, nodes.end() );
10805         }
10806         // set proxy sub-mesh
10807         uvPtVec.resize( nodes.size() );
10808         u2n = nodes.begin();
10809         BRepAdaptor_Curve2d curve( shr->GeomEdge(), F );
10810         for ( size_t i = 0; i < nodes.size(); ++i, ++u2n )
10811         {
10812           uvPtVec[ i ].node = u2n->second;
10813           uvPtVec[ i ].param = u2n->first;
10814           uvPtVec[ i ].SetUV( curve.Value( u2n->first ).XY() );
10815         }
10816         StdMeshers_FaceSide fSide( uvPtVec, F, shr->GeomEdge(), _mesh );
10817         StdMeshers_ViscousLayers2D::SetProxyMeshOfEdge( fSide );
10818       }
10819
10820       // set proxy mesh of EDGEs with layers
10821       vector< _LayerEdge* > edges;
10822       for ( size_t iS = 0; iS < subEOS.size(); ++iS )
10823       {
10824         _EdgesOnShape& eos = * subEOS[ iS ];
10825         if ( eos.ShapeType() != TopAbs_EDGE ) continue;
10826
10827         const TopoDS_Edge& E = TopoDS::Edge( eos._shape );
10828         data.SortOnEdge( E, eos._edges );
10829
10830         edges.clear();
10831         if ( _EdgesOnShape* eov = data.GetShapeEdges( helper.IthVertex( 0, E, /*CumOri=*/false )))
10832           if ( !eov->_edges.empty() )
10833             edges.push_back( eov->_edges[0] ); // on 1st VERTEX
10834
10835         edges.insert( edges.end(), eos._edges.begin(), eos._edges.end() );
10836
10837         if ( _EdgesOnShape* eov = data.GetShapeEdges( helper.IthVertex( 1, E, /*CumOri=*/false )))
10838           if ( !eov->_edges.empty() )
10839             edges.push_back( eov->_edges[0] ); // on last VERTEX
10840
10841         uvPtVec.resize( edges.size() );
10842         for ( size_t i = 0; i < edges.size(); ++i )
10843         {
10844           uvPtVec[ i ].node = edges[i]->_nodes.back();
10845           uvPtVec[ i ].param = helper.GetNodeU( E, edges[i]->_nodes[0] );
10846           uvPtVec[ i ].SetUV( helper.GetNodeUV( F, edges[i]->_nodes.back() ));
10847         }
10848         BRep_Tool::Range( E, uvPtVec[0].param, uvPtVec.back().param );
10849         StdMeshers_FaceSide fSide( uvPtVec, F, E, _mesh );
10850         StdMeshers_ViscousLayers2D::SetProxyMeshOfEdge( fSide );
10851       }
10852       // temporary clear the FACE sub-mesh from faces made by refine()
10853       vector< const SMDS_MeshElement* > elems;
10854       elems.reserve( smDS->NbElements() + smDS->NbNodes() );
10855       for ( SMDS_ElemIteratorPtr ite = smDS->GetElements(); ite->more(); )
10856         elems.push_back( ite->next() );
10857       for ( SMDS_NodeIteratorPtr ite = smDS->GetNodes(); ite->more(); )
10858         elems.push_back( ite->next() );
10859       smDS->Clear();
10860
10861       // compute the mesh on the FACE
10862       sm->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
10863       sm->ComputeStateEngine( SMESH_subMesh::COMPUTE_SUBMESH );
10864
10865       // re-fill proxy sub-meshes of the FACE
10866       for ( size_t i = 0 ; i < _sdVec.size(); ++i )
10867         if (( psm = _sdVec[i]._proxyMesh->getFaceSubM( F )))
10868           for ( SMDS_ElemIteratorPtr ite = smDS->GetElements(); ite->more(); )
10869             psm->AddElement( ite->next() );
10870
10871       // re-fill smDS
10872       for ( size_t i = 0; i < elems.size(); ++i )
10873         smDS->AddElement( elems[i] );
10874
10875       if ( sm->GetComputeState() != SMESH_subMesh::COMPUTE_OK )
10876         return error( errMsg );
10877
10878     } // end of re-meshing in case of failed smoothing
10879     else
10880     {
10881       // No wrongly shaped faces remain; final smooth. Set node XYZ.
10882       bool isStructuredFixed = false;
10883       if ( SMESH_2D_Algo* algo = dynamic_cast<SMESH_2D_Algo*>( sm->GetAlgo() ))
10884         isStructuredFixed = algo->FixInternalNodes( *data._proxyMesh, F );
10885       if ( !isStructuredFixed )
10886       {
10887         if ( isConcaveFace ) // fix narrow faces by swapping diagonals
10888           fixBadFaces( F, helper, /*is2D=*/false, ++shriStep );
10889
10890         for ( int st = 3; st; --st )
10891         {
10892           switch( st ) {
10893           case 1: smoothType = _SmoothNode::LAPLACIAN; break;
10894           case 2: smoothType = _SmoothNode::LAPLACIAN; break;
10895           case 3: smoothType = _SmoothNode::ANGULAR; break;
10896           }
10897           dumpFunction(SMESH_Comment("shrinkFace")<<f2sd->first<<"_st"<<++smooStep); // debug
10898           for ( size_t i = 0; i < nodesToSmooth.size(); ++i )
10899           {
10900             nodesToSmooth[i].Smooth( nbBad,surface,helper,refSign,
10901                                      smoothType,/*set3D=*/st==1 );
10902           }
10903           dumpFunctionEnd();
10904         }
10905       }
10906       if ( !getMeshDS()->IsEmbeddedMode() )
10907         // Log node movement
10908         for ( size_t i = 0; i < nodesToSmooth.size(); ++i )
10909         {
10910           SMESH_TNodeXYZ p ( nodesToSmooth[i]._node );
10911           getMeshDS()->MoveNode( nodesToSmooth[i]._node, p.X(), p.Y(), p.Z() );
10912         }
10913     }
10914
10915     // Set an event listener to clear FACE sub-mesh together with SOLID sub-mesh
10916     VISCOUS_3D::ToClearSubWithMain( sm, data._solid );
10917     if ( data2 )
10918       VISCOUS_3D::ToClearSubWithMain( sm, data2->_solid );
10919
10920   } // loop on FACES to shrink mesh on
10921
10922
10923   // Replace source nodes by target nodes in shrinked mesh edges
10924
10925   map< int, _Shrinker1D >::iterator e2shr = e2shrMap.begin();
10926   for ( ; e2shr != e2shrMap.end(); ++e2shr )
10927     e2shr->second.SwapSrcTgtNodes( getMeshDS() );
10928
10929   return true;
10930 }
10931
10932 //================================================================================
10933 /*!
10934  * \brief Computes 2d shrink direction and finds nodes limiting shrinking
10935  */
10936 //================================================================================
10937
10938 bool _ViscousBuilder::prepareEdgeToShrink( _LayerEdge&            edge,
10939                                            _EdgesOnShape&         eos,
10940                                            SMESH_MesherHelper&    helper,
10941                                            const SMESHDS_SubMesh* faceSubMesh)
10942 {
10943   const SMDS_MeshNode* srcNode = edge._nodes[0];
10944   const SMDS_MeshNode* tgtNode = edge._nodes.back();
10945
10946   if ( eos.SWOLType() == TopAbs_FACE )
10947   {
10948     if ( tgtNode->GetPosition()->GetDim() != 2 ) // not inflated edge
10949     {
10950       edge._pos.clear();
10951       edge.Set( _LayerEdge::SHRUNK );
10952       return srcNode == tgtNode;
10953     }
10954     gp_XY srcUV ( edge._pos[0].X(), edge._pos[0].Y() );          //helper.GetNodeUV( F, srcNode );
10955     gp_XY tgtUV = edge.LastUV( TopoDS::Face( eos._sWOL ), eos ); //helper.GetNodeUV( F, tgtNode );
10956     gp_Vec2d uvDir( srcUV, tgtUV );
10957     double uvLen = uvDir.Magnitude();
10958     uvDir /= uvLen;
10959     edge._normal.SetCoord( uvDir.X(),uvDir.Y(), 0 );
10960     edge._len = uvLen;
10961
10962     //edge._pos.resize(1);
10963     edge._pos[0].SetCoord( tgtUV.X(), tgtUV.Y(), 0 );
10964
10965     // set UV of source node to target node
10966     SMDS_FacePositionPtr pos = tgtNode->GetPosition();
10967     pos->SetUParameter( srcUV.X() );
10968     pos->SetVParameter( srcUV.Y() );
10969   }
10970   else // _sWOL is TopAbs_EDGE
10971   {
10972     if ( tgtNode->GetPosition()->GetDim() != 1 ) // not inflated edge
10973     {
10974       edge._pos.clear();
10975       edge.Set( _LayerEdge::SHRUNK );
10976       return srcNode == tgtNode;
10977     }
10978     const TopoDS_Edge&    E = TopoDS::Edge( eos._sWOL );
10979     SMESHDS_SubMesh* edgeSM = getMeshDS()->MeshElements( E );
10980     if ( !edgeSM || edgeSM->NbElements() == 0 )
10981       return error(SMESH_Comment("Not meshed EDGE ") << getMeshDS()->ShapeToIndex( E ));
10982
10983     const SMDS_MeshNode* n2 = 0;
10984     SMDS_ElemIteratorPtr eIt = srcNode->GetInverseElementIterator(SMDSAbs_Edge);
10985     while ( eIt->more() && !n2 )
10986     {
10987       const SMDS_MeshElement* e = eIt->next();
10988       if ( !edgeSM->Contains(e)) continue;
10989       n2 = e->GetNode( 0 );
10990       if ( n2 == srcNode ) n2 = e->GetNode( 1 );
10991     }
10992     if ( !n2 )
10993       return error(SMESH_Comment("Wrongly meshed EDGE ") << getMeshDS()->ShapeToIndex( E ));
10994
10995     double uSrc = helper.GetNodeU( E, srcNode, n2 );
10996     double uTgt = helper.GetNodeU( E, tgtNode, srcNode );
10997     double u2   = helper.GetNodeU( E, n2,      srcNode );
10998
10999     //edge._pos.clear();
11000
11001     if ( fabs( uSrc-uTgt ) < 0.99 * fabs( uSrc-u2 ))
11002     {
11003       // tgtNode is located so that it does not make faces with wrong orientation
11004       edge.Set( _LayerEdge::SHRUNK );
11005       return true;
11006     }
11007     //edge._pos.resize(1);
11008     edge._pos[0].SetCoord( U_TGT, uTgt );
11009     edge._pos[0].SetCoord( U_SRC, uSrc );
11010     edge._pos[0].SetCoord( LEN_TGT, fabs( uSrc-uTgt ));
11011
11012     edge._simplices.resize( 1 );
11013     edge._simplices[0]._nPrev = n2;
11014
11015     // set U of source node to the target node
11016     SMDS_EdgePositionPtr pos = tgtNode->GetPosition();
11017     pos->SetUParameter( uSrc );
11018   }
11019   return true;
11020 }
11021
11022 //================================================================================
11023 /*!
11024  * \brief Restore position of a sole node of a _LayerEdge based on _noShrinkShapes
11025  */
11026 //================================================================================
11027
11028 void _ViscousBuilder::restoreNoShrink( _LayerEdge& edge ) const
11029 {
11030   if ( edge._nodes.size() == 1 )
11031   {
11032     edge._pos.clear();
11033     edge._len = 0;
11034
11035     const SMDS_MeshNode* srcNode = edge._nodes[0];
11036     TopoDS_Shape S = SMESH_MesherHelper::GetSubShapeByNode( srcNode, getMeshDS() );
11037     if ( S.IsNull() ) return;
11038
11039     gp_Pnt p;
11040
11041     switch ( S.ShapeType() )
11042     {
11043     case TopAbs_EDGE:
11044     {
11045       double f,l;
11046       TopLoc_Location loc;
11047       Handle(Geom_Curve) curve = BRep_Tool::Curve( TopoDS::Edge( S ), loc, f, l );
11048       if ( curve.IsNull() ) return;
11049       SMDS_EdgePositionPtr ePos = srcNode->GetPosition();
11050       p = curve->Value( ePos->GetUParameter() );
11051       break;
11052     }
11053     case TopAbs_VERTEX:
11054     {
11055       p = BRep_Tool::Pnt( TopoDS::Vertex( S ));
11056       break;
11057     }
11058     default: return;
11059     }
11060     getMeshDS()->MoveNode( srcNode, p.X(), p.Y(), p.Z() );
11061     dumpMove( srcNode );
11062   }
11063 }
11064
11065 //================================================================================
11066 /*!
11067  * \brief Try to fix triangles with high aspect ratio by swapping diagonals
11068  */
11069 //================================================================================
11070
11071 void _ViscousBuilder::fixBadFaces(const TopoDS_Face&          F,
11072                                   SMESH_MesherHelper&         helper,
11073                                   const bool                  is2D,
11074                                   const int                   step,
11075                                   set<const SMDS_MeshNode*> * involvedNodes)
11076 {
11077   SMESH::Controls::AspectRatio qualifier;
11078   SMESH::Controls::TSequenceOfXYZ points(3), points1(3), points2(3);
11079   const double maxAspectRatio = is2D ? 4. : 2;
11080   _NodeCoordHelper xyz( F, helper, is2D );
11081
11082   // find bad triangles
11083
11084   vector< const SMDS_MeshElement* > badTrias;
11085   vector< double >                  badAspects;
11086   SMESHDS_SubMesh*      sm = helper.GetMeshDS()->MeshElements( F );
11087   SMDS_ElemIteratorPtr fIt = sm->GetElements();
11088   while ( fIt->more() )
11089   {
11090     const SMDS_MeshElement * f = fIt->next();
11091     if ( f->NbCornerNodes() != 3 ) continue;
11092     for ( int iP = 0; iP < 3; ++iP ) points(iP+1) = xyz( f->GetNode(iP));
11093     double aspect = qualifier.GetValue( points );
11094     if ( aspect > maxAspectRatio )
11095     {
11096       badTrias.push_back( f );
11097       badAspects.push_back( aspect );
11098     }
11099   }
11100   if ( step == 1 )
11101   {
11102     dumpFunction(SMESH_Comment("beforeSwapDiagonals_F")<<helper.GetSubShapeID());
11103     SMDS_ElemIteratorPtr fIt = sm->GetElements();
11104     while ( fIt->more() )
11105     {
11106       const SMDS_MeshElement * f = fIt->next();
11107       if ( f->NbCornerNodes() == 3 )
11108         dumpChangeNodes( f );
11109     }
11110     dumpFunctionEnd();
11111   }
11112   if ( badTrias.empty() )
11113     return;
11114
11115   // find couples of faces to swap diagonal
11116
11117   typedef pair < const SMDS_MeshElement* , const SMDS_MeshElement* > T2Trias;
11118   vector< T2Trias > triaCouples; 
11119
11120   TIDSortedElemSet involvedFaces, emptySet;
11121   for ( size_t iTia = 0; iTia < badTrias.size(); ++iTia )
11122   {
11123     T2Trias trias    [3];
11124     double  aspRatio [3];
11125     int i1, i2, i3;
11126
11127     if ( !involvedFaces.insert( badTrias[iTia] ).second )
11128       continue;
11129     for ( int iP = 0; iP < 3; ++iP )
11130       points(iP+1) = xyz( badTrias[iTia]->GetNode(iP));
11131
11132     // find triangles adjacent to badTrias[iTia] with better aspect ratio after diag-swaping
11133     int bestCouple = -1;
11134     for ( int iSide = 0; iSide < 3; ++iSide )
11135     {
11136       const SMDS_MeshNode* n1 = badTrias[iTia]->GetNode( iSide );
11137       const SMDS_MeshNode* n2 = badTrias[iTia]->GetNode(( iSide+1 ) % 3 );
11138       trias [iSide].first  = badTrias[iTia];
11139       trias [iSide].second = SMESH_MeshAlgos::FindFaceInSet( n1, n2, emptySet, involvedFaces,
11140                                                              & i1, & i2 );
11141       if (( ! trias[iSide].second ) ||
11142           ( trias[iSide].second->NbCornerNodes() != 3 ) ||
11143           ( ! sm->Contains( trias[iSide].second )))
11144         continue;
11145
11146       // aspect ratio of an adjacent tria
11147       for ( int iP = 0; iP < 3; ++iP )
11148         points2(iP+1) = xyz( trias[iSide].second->GetNode(iP));
11149       double aspectInit = qualifier.GetValue( points2 );
11150
11151       // arrange nodes as after diag-swaping
11152       if ( helper.WrapIndex( i1+1, 3 ) == i2 )
11153         i3 = helper.WrapIndex( i1-1, 3 );
11154       else
11155         i3 = helper.WrapIndex( i1+1, 3 );
11156       points1 = points;
11157       points1( 1+ iSide ) = points2( 1+ i3 );
11158       points2( 1+ i2    ) = points1( 1+ ( iSide+2 ) % 3 );
11159
11160       // aspect ratio after diag-swaping
11161       aspRatio[ iSide ] = qualifier.GetValue( points1 ) + qualifier.GetValue( points2 );
11162       if ( aspRatio[ iSide ] > aspectInit + badAspects[ iTia ] )
11163         continue;
11164
11165       // prevent inversion of a triangle
11166       gp_Vec norm1 = gp_Vec( points1(1), points1(3) ) ^ gp_Vec( points1(1), points1(2) );
11167       gp_Vec norm2 = gp_Vec( points2(1), points2(3) ) ^ gp_Vec( points2(1), points2(2) );
11168       if ( norm1 * norm2 < 0. && norm1.Angle( norm2 ) > 70./180.*M_PI )
11169         continue;
11170
11171       if ( bestCouple < 0 || aspRatio[ bestCouple ] > aspRatio[ iSide ] )
11172         bestCouple = iSide;
11173     }
11174
11175     if ( bestCouple >= 0 )
11176     {
11177       triaCouples.push_back( trias[bestCouple] );
11178       involvedFaces.insert ( trias[bestCouple].second );
11179     }
11180     else
11181     {
11182       involvedFaces.erase( badTrias[iTia] );
11183     }
11184   }
11185   if ( triaCouples.empty() )
11186     return;
11187
11188   // swap diagonals
11189
11190   SMESH_MeshEditor editor( helper.GetMesh() );
11191   dumpFunction(SMESH_Comment("beforeSwapDiagonals_F")<<helper.GetSubShapeID()<<"_"<<step);
11192   for ( size_t i = 0; i < triaCouples.size(); ++i )
11193   {
11194     dumpChangeNodes( triaCouples[i].first );
11195     dumpChangeNodes( triaCouples[i].second );
11196     editor.InverseDiag( triaCouples[i].first, triaCouples[i].second );
11197   }
11198
11199   if ( involvedNodes )
11200     for ( size_t i = 0; i < triaCouples.size(); ++i )
11201     {
11202       involvedNodes->insert( triaCouples[i].first->begin_nodes(),
11203                              triaCouples[i].first->end_nodes() );
11204       involvedNodes->insert( triaCouples[i].second->begin_nodes(),
11205                              triaCouples[i].second->end_nodes() );
11206     }
11207
11208   // just for debug dump resulting triangles
11209   dumpFunction(SMESH_Comment("swapDiagonals_F")<<helper.GetSubShapeID()<<"_"<<step);
11210   for ( size_t i = 0; i < triaCouples.size(); ++i )
11211   {
11212     dumpChangeNodes( triaCouples[i].first );
11213     dumpChangeNodes( triaCouples[i].second );
11214   }
11215 }
11216
11217 //================================================================================
11218 /*!
11219  * \brief Move target node to it's final position on the FACE during shrinking
11220  */
11221 //================================================================================
11222
11223 bool _LayerEdge::SetNewLength2d( Handle(Geom_Surface)& surface,
11224                                  const TopoDS_Face&    F,
11225                                  _EdgesOnShape&        eos,
11226                                  SMESH_MesherHelper&   helper )
11227 {
11228   if ( Is( SHRUNK ))
11229     return false; // already at the target position
11230
11231   SMDS_MeshNode* tgtNode = const_cast< SMDS_MeshNode*& >( _nodes.back() );
11232
11233   if ( eos.SWOLType() == TopAbs_FACE )
11234   {
11235     gp_XY    curUV = helper.GetNodeUV( F, tgtNode );
11236     gp_Pnt2d tgtUV( _pos[0].X(), _pos[0].Y() );
11237     gp_Vec2d uvDir( _normal.X(), _normal.Y() );
11238     const double uvLen = tgtUV.Distance( curUV );
11239     const double kSafe = Max( 0.5, 1. - 0.1 * _simplices.size() );
11240
11241     // Select shrinking step such that not to make faces with wrong orientation.
11242     double stepSize = 1e100;
11243     for ( size_t i = 0; i < _simplices.size(); ++i )
11244     {
11245       if ( !_simplices[i]._nPrev->isMarked() ||
11246            !_simplices[i]._nNext->isMarked() )
11247         continue; // simplex of quadrangle created by addBoundaryElements()
11248
11249       // find intersection of 2 lines: curUV-tgtUV and that connecting simplex nodes
11250       gp_XY uvN1 = helper.GetNodeUV( F, _simplices[i]._nPrev );
11251       gp_XY uvN2 = helper.GetNodeUV( F, _simplices[i]._nNext );
11252       gp_XY dirN = uvN2 - uvN1;
11253       double det = uvDir.Crossed( dirN );
11254       if ( Abs( det )  < std::numeric_limits<double>::min() ) continue;
11255       gp_XY dirN2Cur = curUV - uvN1;
11256       double step = dirN.Crossed( dirN2Cur ) / det;
11257       if ( step > 0 )
11258         stepSize = Min( step, stepSize );
11259     }
11260     gp_Pnt2d newUV;
11261     if ( uvLen <= stepSize )
11262     {
11263       newUV = tgtUV;
11264       Set( SHRUNK );
11265       //_pos.clear();
11266     }
11267     else if ( stepSize > 0 )
11268     {
11269       newUV = curUV + uvDir.XY() * stepSize * kSafe;
11270     }
11271     else
11272     {
11273       return true;
11274     }
11275     SMDS_FacePositionPtr pos = tgtNode->GetPosition();
11276     pos->SetUParameter( newUV.X() );
11277     pos->SetVParameter( newUV.Y() );
11278
11279 #ifdef __myDEBUG
11280     gp_Pnt p = surface->Value( newUV.X(), newUV.Y() );
11281     tgtNode->setXYZ( p.X(), p.Y(), p.Z() );
11282     dumpMove( tgtNode );
11283 #endif
11284   }
11285   else // _sWOL is TopAbs_EDGE
11286   {
11287     const TopoDS_Edge&      E = TopoDS::Edge( eos._sWOL );
11288     const SMDS_MeshNode*   n2 = _simplices[0]._nPrev;
11289     SMDS_EdgePositionPtr tgtPos = tgtNode->GetPosition();
11290
11291     const double u2     = helper.GetNodeU( E, n2, tgtNode );
11292     const double uSrc   = _pos[0].Coord( U_SRC );
11293     const double lenTgt = _pos[0].Coord( LEN_TGT );
11294
11295     double newU = _pos[0].Coord( U_TGT );
11296     if ( lenTgt < 0.99 * fabs( uSrc-u2 )) // n2 got out of src-tgt range
11297     {
11298       Set( _LayerEdge::SHRUNK );
11299       //_pos.clear();
11300     }
11301     else
11302     {
11303       newU = 0.1 * tgtPos->GetUParameter() + 0.9 * u2;
11304     }
11305     tgtPos->SetUParameter( newU );
11306 #ifdef __myDEBUG
11307     gp_XY newUV = helper.GetNodeUV( F, tgtNode, _nodes[0]);
11308     gp_Pnt p = surface->Value( newUV.X(), newUV.Y() );
11309     tgtNode->setXYZ( p.X(), p.Y(), p.Z() );
11310     dumpMove( tgtNode );
11311 #endif
11312   }
11313
11314   return true;
11315 }
11316
11317 //================================================================================
11318 /*!
11319  * \brief Perform smooth on the FACE
11320  *  \retval bool - true if the node has been moved
11321  */
11322 //================================================================================
11323
11324 bool _SmoothNode::Smooth(int&                  nbBad,
11325                          Handle(Geom_Surface)& surface,
11326                          SMESH_MesherHelper&   helper,
11327                          const double          refSign,
11328                          SmoothType            how,
11329                          bool                  set3D)
11330 {
11331   const TopoDS_Face& face = TopoDS::Face( helper.GetSubShape() );
11332
11333   // get uv of surrounding nodes
11334   vector<gp_XY> uv( _simplices.size() );
11335   for ( size_t i = 0; i < _simplices.size(); ++i )
11336     uv[i] = helper.GetNodeUV( face, _simplices[i]._nPrev, _node );
11337
11338   // compute new UV for the node
11339   gp_XY newPos (0,0);
11340   if ( how == TFI && _simplices.size() == 4 )
11341   {
11342     gp_XY corners[4];
11343     for ( size_t i = 0; i < _simplices.size(); ++i )
11344       if ( _simplices[i]._nOpp )
11345         corners[i] = helper.GetNodeUV( face, _simplices[i]._nOpp, _node );
11346       else
11347         throw SALOME_Exception(LOCALIZED("TFI smoothing: _Simplex::_nOpp not set!"));
11348
11349     newPos = helper.calcTFI ( 0.5, 0.5,
11350                               corners[0], corners[1], corners[2], corners[3],
11351                               uv[1], uv[2], uv[3], uv[0] );
11352   }
11353   else if ( how == ANGULAR )
11354   {
11355     newPos = computeAngularPos( uv, helper.GetNodeUV( face, _node ), refSign );
11356   }
11357   else if ( how == CENTROIDAL && _simplices.size() > 3 )
11358   {
11359     // average centers of diagonals wieghted with their reciprocal lengths
11360     if ( _simplices.size() == 4 )
11361     {
11362       double w1 = 1. / ( uv[2]-uv[0] ).SquareModulus();
11363       double w2 = 1. / ( uv[3]-uv[1] ).SquareModulus();
11364       newPos = ( w1 * ( uv[2]+uv[0] ) + w2 * ( uv[3]+uv[1] )) / ( w1+w2 ) / 2;
11365     }
11366     else
11367     {
11368       double sumWeight = 0;
11369       int nb = _simplices.size() == 4 ? 2 : _simplices.size();
11370       for ( int i = 0; i < nb; ++i )
11371       {
11372         int iFrom = i + 2;
11373         int iTo   = i + _simplices.size() - 1;
11374         for ( int j = iFrom; j < iTo; ++j )
11375         {
11376           int i2 = SMESH_MesherHelper::WrapIndex( j, _simplices.size() );
11377           double w = 1. / ( uv[i]-uv[i2] ).SquareModulus();
11378           sumWeight += w;
11379           newPos += w * ( uv[i]+uv[i2] );
11380         }
11381       }
11382       newPos /= 2 * sumWeight; // 2 is to get a middle between uv's
11383     }
11384   }
11385   else
11386   {
11387     // Laplacian smooth
11388     for ( size_t i = 0; i < _simplices.size(); ++i )
11389       newPos += uv[i];
11390     newPos /= _simplices.size();
11391   }
11392
11393   // count quality metrics (orientation) of triangles around the node
11394   int nbOkBefore = 0;
11395   gp_XY tgtUV = helper.GetNodeUV( face, _node );
11396   for ( size_t i = 0; i < _simplices.size(); ++i )
11397     nbOkBefore += _simplices[i].IsForward( tgtUV, _node, face, helper, refSign );
11398
11399   int nbOkAfter = 0;
11400   for ( size_t i = 0; i < _simplices.size(); ++i )
11401     nbOkAfter += _simplices[i].IsForward( newPos, _node, face, helper, refSign );
11402
11403   if ( nbOkAfter < nbOkBefore )
11404   {
11405     nbBad += _simplices.size() - nbOkBefore;
11406     return false;
11407   }
11408
11409   SMDS_FacePositionPtr pos = _node->GetPosition();
11410   pos->SetUParameter( newPos.X() );
11411   pos->SetVParameter( newPos.Y() );
11412
11413 #ifdef __myDEBUG
11414   set3D = true;
11415 #endif
11416   if ( set3D )
11417   {
11418     gp_Pnt p = surface->Value( newPos.X(), newPos.Y() );
11419     const_cast< SMDS_MeshNode* >( _node )->setXYZ( p.X(), p.Y(), p.Z() );
11420     dumpMove( _node );
11421   }
11422
11423   nbBad += _simplices.size() - nbOkAfter;
11424   return ( (tgtUV-newPos).SquareModulus() > 1e-10 );
11425 }
11426
11427 //================================================================================
11428 /*!
11429  * \brief Computes new UV using angle based smoothing technique
11430  */
11431 //================================================================================
11432
11433 gp_XY _SmoothNode::computeAngularPos(vector<gp_XY>& uv,
11434                                      const gp_XY&   uvToFix,
11435                                      const double   refSign)
11436 {
11437   uv.push_back( uv.front() );
11438
11439   vector< gp_XY >  edgeDir ( uv.size() );
11440   vector< double > edgeSize( uv.size() );
11441   for ( size_t i = 1; i < edgeDir.size(); ++i )
11442   {
11443     edgeDir [i-1] = uv[i] - uv[i-1];
11444     edgeSize[i-1] = edgeDir[i-1].Modulus();
11445     if ( edgeSize[i-1] < numeric_limits<double>::min() )
11446       edgeDir[i-1].SetX( 100 );
11447     else
11448       edgeDir[i-1] /= edgeSize[i-1] * refSign;
11449   }
11450   edgeDir.back()  = edgeDir.front();
11451   edgeSize.back() = edgeSize.front();
11452
11453   gp_XY  newPos(0,0);
11454   //int    nbEdges = 0;
11455   double sumSize = 0;
11456   for ( size_t i = 1; i < edgeDir.size(); ++i )
11457   {
11458     if ( edgeDir[i-1].X() > 1. ) continue;
11459     int i1 = i-1;
11460     while ( edgeDir[i].X() > 1. && ++i < edgeDir.size() );
11461     if ( i == edgeDir.size() ) break;
11462     gp_XY p = uv[i];
11463     gp_XY norm1( -edgeDir[i1].Y(), edgeDir[i1].X() );
11464     gp_XY norm2( -edgeDir[i].Y(),  edgeDir[i].X() );
11465     gp_XY bisec = norm1 + norm2;
11466     double bisecSize = bisec.Modulus();
11467     if ( bisecSize < numeric_limits<double>::min() )
11468     {
11469       bisec = -edgeDir[i1] + edgeDir[i];
11470       bisecSize = bisec.Modulus();
11471     }
11472     bisec /= bisecSize;
11473
11474     gp_XY  dirToN  = uvToFix - p;
11475     double distToN = dirToN.Modulus();
11476     if ( bisec * dirToN < 0 )
11477       distToN = -distToN;
11478
11479     newPos += ( p + bisec * distToN ) * ( edgeSize[i1] + edgeSize[i] );
11480     //++nbEdges;
11481     sumSize += edgeSize[i1] + edgeSize[i];
11482   }
11483   newPos /= /*nbEdges * */sumSize;
11484   return newPos;
11485 }
11486
11487 //================================================================================
11488 /*!
11489  * \brief Delete _SolidData
11490  */
11491 //================================================================================
11492
11493 _SolidData::~_SolidData()
11494 {
11495   TNode2Edge::iterator n2e = _n2eMap.begin();
11496   for ( ; n2e != _n2eMap.end(); ++n2e )
11497   {
11498     _LayerEdge* & e = n2e->second;
11499     if ( e )
11500     {
11501       delete e->_curvature;
11502       if ( e->_2neibors )
11503         delete e->_2neibors->_plnNorm;
11504       delete e->_2neibors;
11505     }
11506     delete e;
11507     e = 0;
11508   }
11509   _n2eMap.clear();
11510
11511   delete _helper;
11512   _helper = 0;
11513 }
11514
11515 //================================================================================
11516 /*!
11517  * \brief Keep a _LayerEdge inflated along the EDGE
11518  */
11519 //================================================================================
11520
11521 void _Shrinker1D::AddEdge( const _LayerEdge*   e,
11522                            _EdgesOnShape&      eos,
11523                            SMESH_MesherHelper& helper )
11524 {
11525   // init
11526   if ( _nodes.empty() )
11527   {
11528     _edges[0] = _edges[1] = 0;
11529     _done = false;
11530   }
11531   // check _LayerEdge
11532   if ( e == _edges[0] || e == _edges[1] || e->_nodes.size() < 2 )
11533     return;
11534   if ( eos.SWOLType() != TopAbs_EDGE )
11535     throw SALOME_Exception(LOCALIZED("Wrong _LayerEdge is added"));
11536   if ( _edges[0] && !_geomEdge.IsSame( eos._sWOL ))
11537     throw SALOME_Exception(LOCALIZED("Wrong _LayerEdge is added"));
11538
11539   // store _LayerEdge
11540   _geomEdge = TopoDS::Edge( eos._sWOL );
11541   double f,l;
11542   BRep_Tool::Range( _geomEdge, f,l );
11543   double u = helper.GetNodeU( _geomEdge, e->_nodes[0], e->_nodes.back());
11544   _edges[ u < 0.5*(f+l) ? 0 : 1 ] = e;
11545
11546   // Update _nodes
11547
11548   const SMDS_MeshNode* tgtNode0 = TgtNode( 0 );
11549   const SMDS_MeshNode* tgtNode1 = TgtNode( 1 );
11550
11551   if ( _nodes.empty() )
11552   {
11553     SMESHDS_SubMesh * eSubMesh = helper.GetMeshDS()->MeshElements( _geomEdge );
11554     if ( !eSubMesh || eSubMesh->NbNodes() < 1 )
11555       return;
11556     TopLoc_Location loc;
11557     Handle(Geom_Curve) C = BRep_Tool::Curve( _geomEdge, loc, f,l );
11558     GeomAdaptor_Curve aCurve(C, f,l);
11559     const double totLen = GCPnts_AbscissaPoint::Length(aCurve, f, l);
11560
11561     int nbExpectNodes = eSubMesh->NbNodes();
11562     _initU  .reserve( nbExpectNodes );
11563     _normPar.reserve( nbExpectNodes );
11564     _nodes  .reserve( nbExpectNodes );
11565     SMDS_NodeIteratorPtr nIt = eSubMesh->GetNodes();
11566     while ( nIt->more() )
11567     {
11568       const SMDS_MeshNode* node = nIt->next();
11569
11570       // skip refinement nodes
11571       if ( node->NbInverseElements(SMDSAbs_Edge) == 0 ||
11572            node == tgtNode0 || node == tgtNode1 )
11573         continue;
11574       bool hasMarkedFace = false;
11575       SMDS_ElemIteratorPtr fIt = node->GetInverseElementIterator(SMDSAbs_Face);
11576       while ( fIt->more() && !hasMarkedFace )
11577         hasMarkedFace = fIt->next()->isMarked();
11578       if ( !hasMarkedFace )
11579         continue;
11580
11581       _nodes.push_back( node );
11582       _initU.push_back( helper.GetNodeU( _geomEdge, node ));
11583       double len = GCPnts_AbscissaPoint::Length(aCurve, f, _initU.back());
11584       _normPar.push_back(  len / totLen );
11585     }
11586   }
11587   else
11588   {
11589     // remove target node of the _LayerEdge from _nodes
11590     size_t nbFound = 0;
11591     for ( size_t i = 0; i < _nodes.size(); ++i )
11592       if ( !_nodes[i] || _nodes[i] == tgtNode0 || _nodes[i] == tgtNode1 )
11593         _nodes[i] = 0, nbFound++;
11594     if ( nbFound == _nodes.size() )
11595       _nodes.clear();
11596   }
11597 }
11598
11599 //================================================================================
11600 /*!
11601  * \brief Move nodes on EDGE from ends where _LayerEdge's are inflated
11602  */
11603 //================================================================================
11604
11605 void _Shrinker1D::Compute(bool set3D, SMESH_MesherHelper& helper)
11606 {
11607   if ( _done || _nodes.empty())
11608     return;
11609   const _LayerEdge* e = _edges[0];
11610   if ( !e ) e = _edges[1];
11611   if ( !e ) return;
11612
11613   _done =  (( !_edges[0] || _edges[0]->Is( _LayerEdge::SHRUNK )) &&
11614             ( !_edges[1] || _edges[1]->Is( _LayerEdge::SHRUNK )));
11615
11616   double f,l;
11617   if ( set3D || _done )
11618   {
11619     Handle(Geom_Curve) C = BRep_Tool::Curve(_geomEdge, f,l);
11620     GeomAdaptor_Curve aCurve(C, f,l);
11621
11622     if ( _edges[0] )
11623       f = helper.GetNodeU( _geomEdge, _edges[0]->_nodes.back(), _nodes[0] );
11624     if ( _edges[1] )
11625       l = helper.GetNodeU( _geomEdge, _edges[1]->_nodes.back(), _nodes.back() );
11626     double totLen = GCPnts_AbscissaPoint::Length( aCurve, f, l );
11627
11628     for ( size_t i = 0; i < _nodes.size(); ++i )
11629     {
11630       if ( !_nodes[i] ) continue;
11631       double len = totLen * _normPar[i];
11632       GCPnts_AbscissaPoint discret( aCurve, len, f );
11633       if ( !discret.IsDone() )
11634         return throw SALOME_Exception(LOCALIZED("GCPnts_AbscissaPoint failed"));
11635       double u = discret.Parameter();
11636       SMDS_EdgePositionPtr pos = _nodes[i]->GetPosition();
11637       pos->SetUParameter( u );
11638       gp_Pnt p = C->Value( u );
11639       const_cast< SMDS_MeshNode*>( _nodes[i] )->setXYZ( p.X(), p.Y(), p.Z() );
11640     }
11641   }
11642   else
11643   {
11644     BRep_Tool::Range( _geomEdge, f,l );
11645     if ( _edges[0] )
11646       f = helper.GetNodeU( _geomEdge, _edges[0]->_nodes.back(), _nodes[0] );
11647     if ( _edges[1] )
11648       l = helper.GetNodeU( _geomEdge, _edges[1]->_nodes.back(), _nodes.back() );
11649     
11650     for ( size_t i = 0; i < _nodes.size(); ++i )
11651     {
11652       if ( !_nodes[i] ) continue;
11653       double u = f * ( 1-_normPar[i] ) + l * _normPar[i];
11654       SMDS_EdgePositionPtr pos = _nodes[i]->GetPosition();
11655       pos->SetUParameter( u );
11656     }
11657   }
11658 }
11659
11660 //================================================================================
11661 /*!
11662  * \brief Restore initial parameters of nodes on EDGE
11663  */
11664 //================================================================================
11665
11666 void _Shrinker1D::RestoreParams()
11667 {
11668   if ( _done )
11669     for ( size_t i = 0; i < _nodes.size(); ++i )
11670     {
11671       if ( !_nodes[i] ) continue;
11672       SMDS_EdgePositionPtr pos = _nodes[i]->GetPosition();
11673       pos->SetUParameter( _initU[i] );
11674     }
11675   _done = false;
11676 }
11677
11678 //================================================================================
11679 /*!
11680  * \brief Replace source nodes by target nodes in shrinked mesh edges
11681  */
11682 //================================================================================
11683
11684 void _Shrinker1D::SwapSrcTgtNodes( SMESHDS_Mesh* mesh )
11685 {
11686   const SMDS_MeshNode* nodes[3];
11687   for ( int i = 0; i < 2; ++i )
11688   {
11689     if ( !_edges[i] ) continue;
11690
11691     SMESHDS_SubMesh * eSubMesh = mesh->MeshElements( _geomEdge );
11692     if ( !eSubMesh ) return;
11693     const SMDS_MeshNode* srcNode = _edges[i]->_nodes[0];
11694     const SMDS_MeshNode* tgtNode = _edges[i]->_nodes.back();
11695     const SMDS_MeshNode* scdNode = _edges[i]->_nodes[1];
11696     SMDS_ElemIteratorPtr eIt = srcNode->GetInverseElementIterator(SMDSAbs_Edge);
11697     while ( eIt->more() )
11698     {
11699       const SMDS_MeshElement* e = eIt->next();
11700       if ( !eSubMesh->Contains( e ) || e->GetNodeIndex( scdNode ) >= 0 )
11701           continue;
11702       SMDS_ElemIteratorPtr nIt = e->nodesIterator();
11703       for ( int iN = 0; iN < e->NbNodes(); ++iN )
11704       {
11705         const SMDS_MeshNode* n = static_cast<const SMDS_MeshNode*>( nIt->next() );
11706         nodes[iN] = ( n == srcNode ? tgtNode : n );
11707       }
11708       mesh->ChangeElementNodes( e, nodes, e->NbNodes() );
11709     }
11710   }
11711 }
11712
11713 //================================================================================
11714 /*!
11715  * \brief Creates 2D and 1D elements on boundaries of new prisms
11716  */
11717 //================================================================================
11718
11719 bool _ViscousBuilder::addBoundaryElements(_SolidData& data)
11720 {
11721   SMESH_MesherHelper helper( *_mesh );
11722
11723   vector< const SMDS_MeshNode* > faceNodes;
11724
11725   //for ( size_t i = 0; i < _sdVec.size(); ++i )
11726   {
11727     //_SolidData& data = _sdVec[i];
11728     TopTools_IndexedMapOfShape geomEdges;
11729     TopExp::MapShapes( data._solid, TopAbs_EDGE, geomEdges );
11730     for ( int iE = 1; iE <= geomEdges.Extent(); ++iE )
11731     {
11732       const TopoDS_Edge& E = TopoDS::Edge( geomEdges(iE));
11733       const TGeomID edgeID = getMeshDS()->ShapeToIndex( E );
11734       if ( data._noShrinkShapes.count( edgeID ))
11735         continue;
11736
11737       // Get _LayerEdge's based on E
11738
11739       map< double, const SMDS_MeshNode* > u2nodes;
11740       if ( !SMESH_Algo::GetSortedNodesOnEdge( getMeshDS(), E, /*ignoreMedium=*/false, u2nodes))
11741         continue;
11742
11743       vector< _LayerEdge* > ledges; ledges.reserve( u2nodes.size() );
11744       TNode2Edge & n2eMap = data._n2eMap;
11745       map< double, const SMDS_MeshNode* >::iterator u2n = u2nodes.begin();
11746       {
11747         //check if 2D elements are needed on E
11748         TNode2Edge::iterator n2e = n2eMap.find( u2n->second );
11749         if ( n2e == n2eMap.end() ) continue; // no layers on vertex
11750         ledges.push_back( n2e->second );
11751         u2n++;
11752         if (( n2e = n2eMap.find( u2n->second )) == n2eMap.end() )
11753           continue; // no layers on E
11754         ledges.push_back( n2eMap[ u2n->second ]);
11755
11756         const SMDS_MeshNode* tgtN0 = ledges[0]->_nodes.back();
11757         const SMDS_MeshNode* tgtN1 = ledges[1]->_nodes.back();
11758         int nbSharedPyram = 0;
11759         SMDS_ElemIteratorPtr vIt = tgtN1->GetInverseElementIterator(SMDSAbs_Volume);
11760         while ( vIt->more() )
11761         {
11762           const SMDS_MeshElement* v = vIt->next();
11763           nbSharedPyram += int( v->GetNodeIndex( tgtN0 ) >= 0 );
11764         }
11765         if ( nbSharedPyram > 1 )
11766           continue; // not free border of the pyramid
11767
11768         faceNodes.clear();
11769         faceNodes.push_back( ledges[0]->_nodes[0] );
11770         faceNodes.push_back( ledges[1]->_nodes[0] );
11771         if ( ledges[0]->_nodes.size() > 1 ) faceNodes.push_back( ledges[0]->_nodes[1] );
11772         if ( ledges[1]->_nodes.size() > 1 ) faceNodes.push_back( ledges[1]->_nodes[1] );
11773
11774         if ( getMeshDS()->FindElement( faceNodes, SMDSAbs_Face, /*noMedium=*/true))
11775           continue; // faces already created
11776       }
11777       for ( ++u2n; u2n != u2nodes.end(); ++u2n )
11778         ledges.push_back( n2eMap[ u2n->second ]);
11779
11780       // Find out orientation and type of face to create
11781
11782       bool reverse = false, isOnFace;
11783       TopoDS_Shape F;
11784
11785       map< TGeomID, TopoDS_Shape >::iterator e2f = data._shrinkShape2Shape.find( edgeID );
11786       if (( isOnFace = ( e2f != data._shrinkShape2Shape.end() )))
11787       {
11788         F = e2f->second.Oriented( TopAbs_FORWARD );
11789         reverse = ( helper.GetSubShapeOri( F, E ) == TopAbs_REVERSED );
11790         if ( helper.GetSubShapeOri( data._solid, F ) == TopAbs_REVERSED )
11791           reverse = !reverse, F.Reverse();
11792         if ( helper.IsReversedSubMesh( TopoDS::Face(F) ))
11793           reverse = !reverse;
11794       }
11795       else if ( !data._ignoreFaceIds.count( e2f->first ))
11796       {
11797         // find FACE with layers sharing E
11798         PShapeIteratorPtr fIt = helper.GetAncestors( E, *_mesh, TopAbs_FACE, &data._solid );
11799         if ( fIt->more() )
11800           F = *( fIt->next() );
11801       }
11802       // Find the sub-mesh to add new faces
11803       SMESHDS_SubMesh* sm = 0;
11804       if ( isOnFace )
11805         sm = getMeshDS()->MeshElements( F );
11806       else
11807         sm = data._proxyMesh->getFaceSubM( TopoDS::Face(F), /*create=*/true );
11808       if ( !sm )
11809         return error("error in addBoundaryElements()", data._index);
11810
11811       // Find a proxy sub-mesh of the FACE of an adjacent SOLID, which will use the new boundary
11812       // faces for 3D meshing (PAL23414)
11813       SMESHDS_SubMesh* adjSM = 0;
11814       if ( isOnFace )
11815       {
11816         const TGeomID   faceID = sm->GetID();
11817         PShapeIteratorPtr soIt = helper.GetAncestors( F, *_mesh, TopAbs_SOLID );
11818         while ( const TopoDS_Shape* solid = soIt->next() )
11819           if ( !solid->IsSame( data._solid ))
11820           {
11821             size_t iData = _solids.FindIndex( *solid ) - 1;
11822             if ( iData < _sdVec.size() &&
11823                  _sdVec[ iData ]._ignoreFaceIds.count( faceID ) &&
11824                  _sdVec[ iData ]._shrinkShape2Shape.count( edgeID ) == 0 )
11825             {
11826               SMESH_ProxyMesh::SubMesh* proxySub =
11827                 _sdVec[ iData ]._proxyMesh->getFaceSubM( TopoDS::Face( F ), /*create=*/false);
11828               if ( proxySub && proxySub->NbElements() > 0 )
11829                 adjSM = proxySub;
11830             }
11831           }
11832       }
11833
11834       // Make faces
11835       const int dj1 = reverse ? 0 : 1;
11836       const int dj2 = reverse ? 1 : 0;
11837       vector< const SMDS_MeshElement*> ff; // new faces row
11838       SMESHDS_Mesh* m = getMeshDS();
11839       for ( size_t j = 1; j < ledges.size(); ++j )
11840       {
11841         vector< const SMDS_MeshNode*>&  nn1 = ledges[j-dj1]->_nodes;
11842         vector< const SMDS_MeshNode*>&  nn2 = ledges[j-dj2]->_nodes;
11843         ff.resize( std::max( nn1.size(), nn2.size() ), NULL );
11844         if ( nn1.size() == nn2.size() )
11845         {
11846           if ( isOnFace )
11847             for ( size_t z = 1; z < nn1.size(); ++z )
11848               sm->AddElement( ff[z-1] = m->AddFace( nn1[z-1], nn2[z-1], nn2[z], nn1[z] ));
11849           else
11850             for ( size_t z = 1; z < nn1.size(); ++z )
11851               sm->AddElement( new SMDS_FaceOfNodes( nn1[z-1], nn2[z-1], nn2[z], nn1[z] ));
11852         }
11853         else if ( nn1.size() == 1 )
11854         {
11855           if ( isOnFace )
11856             for ( size_t z = 1; z < nn2.size(); ++z )
11857               sm->AddElement( ff[z-1] = m->AddFace( nn1[0], nn2[z-1], nn2[z] ));
11858           else
11859             for ( size_t z = 1; z < nn2.size(); ++z )
11860               sm->AddElement( new SMDS_FaceOfNodes( nn1[0], nn2[z-1], nn2[z] ));
11861         }
11862         else
11863         {
11864           if ( isOnFace )
11865             for ( size_t z = 1; z < nn1.size(); ++z )
11866               sm->AddElement( ff[z-1] = m->AddFace( nn1[z-1], nn2[0], nn1[z] ));
11867           else
11868             for ( size_t z = 1; z < nn1.size(); ++z )
11869               sm->AddElement( new SMDS_FaceOfNodes( nn1[z-1], nn2[0], nn2[z] ));
11870         }
11871
11872         if ( adjSM ) // add faces to a proxy SM of the adjacent SOLID
11873         {
11874           for ( size_t z = 0; z < ff.size(); ++z )
11875             if ( ff[ z ])
11876               adjSM->AddElement( ff[ z ]);
11877           ff.clear();
11878         }
11879       }
11880
11881       // Make edges
11882       for ( int isFirst = 0; isFirst < 2; ++isFirst )
11883       {
11884         _LayerEdge* edge = isFirst ? ledges.front() : ledges.back();
11885         _EdgesOnShape* eos = data.GetShapeEdges( edge );
11886         if ( eos && eos->SWOLType() == TopAbs_EDGE )
11887         {
11888           vector< const SMDS_MeshNode*>&  nn = edge->_nodes;
11889           if ( nn.size() < 2 || nn[1]->NbInverseElements( SMDSAbs_Edge ) >= 2 )
11890             continue;
11891           helper.SetSubShape( eos->_sWOL );
11892           helper.SetElementsOnShape( true );
11893           for ( size_t z = 1; z < nn.size(); ++z )
11894             helper.AddEdge( nn[z-1], nn[z] );
11895         }
11896       }
11897
11898     } // loop on EDGE's
11899   } // loop on _SolidData's
11900
11901   return true;
11902 }