Salome HOME
0c40b47438c6ccf930cbf9e6aed0dfbfdd3ae30b
[modules/smesh.git] / src / StdMeshers / StdMeshers_ViscousLayers.cxx
1 // Copyright (C) 2007-2016  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 // File      : StdMeshers_ViscousLayers.cxx
21 // Created   : Wed Dec  1 15:15:34 2010
22 // Author    : Edward AGAPOV (eap)
23
24 #include "StdMeshers_ViscousLayers.hxx"
25
26 #include "SMDS_EdgePosition.hxx"
27 #include "SMDS_FaceOfNodes.hxx"
28 #include "SMDS_FacePosition.hxx"
29 #include "SMDS_MeshNode.hxx"
30 #include "SMDS_SetIterator.hxx"
31 #include "SMESHDS_Group.hxx"
32 #include "SMESHDS_Hypothesis.hxx"
33 #include "SMESH_Algo.hxx"
34 #include "SMESH_ComputeError.hxx"
35 #include "SMESH_ControlsDef.hxx"
36 #include "SMESH_Gen.hxx"
37 #include "SMESH_Group.hxx"
38 #include "SMESH_HypoFilter.hxx"
39 #include "SMESH_Mesh.hxx"
40 #include "SMESH_MeshAlgos.hxx"
41 #include "SMESH_MesherHelper.hxx"
42 #include "SMESH_ProxyMesh.hxx"
43 #include "SMESH_subMesh.hxx"
44 #include "SMESH_subMeshEventListener.hxx"
45 #include "StdMeshers_FaceSide.hxx"
46 #include "StdMeshers_ViscousLayers2D.hxx"
47
48 #include <Adaptor3d_HSurface.hxx>
49 #include <BRepAdaptor_Curve.hxx>
50 #include <BRepAdaptor_Curve2d.hxx>
51 #include <BRepAdaptor_Surface.hxx>
52 #include <BRepLProp_SLProps.hxx>
53 #include <BRepOffsetAPI_MakeOffsetShape.hxx>
54 #include <BRep_Tool.hxx>
55 #include <Bnd_B2d.hxx>
56 #include <Bnd_B3d.hxx>
57 #include <ElCLib.hxx>
58 #include <GCPnts_AbscissaPoint.hxx>
59 #include <GCPnts_TangentialDeflection.hxx>
60 #include <Geom2d_Circle.hxx>
61 #include <Geom2d_Line.hxx>
62 #include <Geom2d_TrimmedCurve.hxx>
63 #include <GeomAdaptor_Curve.hxx>
64 #include <GeomLib.hxx>
65 #include <Geom_Circle.hxx>
66 #include <Geom_Curve.hxx>
67 #include <Geom_Line.hxx>
68 #include <Geom_TrimmedCurve.hxx>
69 #include <Precision.hxx>
70 #include <Standard_ErrorHandler.hxx>
71 #include <Standard_Failure.hxx>
72 #include <TColStd_Array1OfReal.hxx>
73 #include <TopExp.hxx>
74 #include <TopExp_Explorer.hxx>
75 #include <TopTools_IndexedMapOfShape.hxx>
76 #include <TopTools_ListOfShape.hxx>
77 #include <TopTools_MapIteratorOfMapOfShape.hxx>
78 #include <TopTools_MapOfShape.hxx>
79 #include <TopoDS.hxx>
80 #include <TopoDS_Edge.hxx>
81 #include <TopoDS_Face.hxx>
82 #include <TopoDS_Vertex.hxx>
83 #include <gp_Ax1.hxx>
84 #include <gp_Cone.hxx>
85 #include <gp_Sphere.hxx>
86 #include <gp_Vec.hxx>
87 #include <gp_XY.hxx>
88
89 #include <cmath>
90 #include <limits>
91 #include <list>
92 #include <queue>
93 #include <string>
94
95 #ifdef _DEBUG_
96 #define __myDEBUG
97 //#define __NOT_INVALIDATE_BAD_SMOOTH
98 #endif
99
100 #define INCREMENTAL_SMOOTH // smooth only if min angle is too small
101 #define BLOCK_INFLATION // of individual _LayerEdge's
102 #define OLD_NEF_POLYGON
103
104 using namespace std;
105
106 //================================================================================
107 namespace VISCOUS_3D
108 {
109   typedef int TGeomID;
110
111   enum UIndex { U_TGT = 1, U_SRC, LEN_TGT };
112
113   const double theMinSmoothCosin = 0.1;
114   const double theSmoothThickToElemSizeRatio = 0.3;
115   const double theMinSmoothTriaAngle = 30;
116   const double theMinSmoothQuadAngle = 45;
117
118   // what part of thickness is allowed till intersection
119   // (defined by SALOME_TESTS/Grids/smesh/viscous_layers_00/A5)
120   const double theThickToIntersection = 1.5;
121
122   bool needSmoothing( double cosin, double tgtThick, double elemSize )
123   {
124     return cosin * tgtThick > theSmoothThickToElemSizeRatio * elemSize;
125   }
126   double getSmoothingThickness( double cosin, double elemSize )
127   {
128     return theSmoothThickToElemSizeRatio * elemSize / cosin;
129   }
130
131   /*!
132    * \brief SMESH_ProxyMesh computed by _ViscousBuilder for a SOLID.
133    * It is stored in a SMESH_subMesh of the SOLID as SMESH_subMeshEventListenerData
134    */
135   struct _MeshOfSolid : public SMESH_ProxyMesh,
136                         public SMESH_subMeshEventListenerData
137   {
138     bool                  _n2nMapComputed;
139     SMESH_ComputeErrorPtr _warning;
140
141     _MeshOfSolid( SMESH_Mesh* mesh)
142       :SMESH_subMeshEventListenerData( /*isDeletable=*/true),_n2nMapComputed(false)
143     {
144       SMESH_ProxyMesh::setMesh( *mesh );
145     }
146
147     // returns submesh for a geom face
148     SMESH_ProxyMesh::SubMesh* getFaceSubM(const TopoDS_Face& F, bool create=false)
149     {
150       TGeomID i = SMESH_ProxyMesh::shapeIndex(F);
151       return create ? SMESH_ProxyMesh::getProxySubMesh(i) : findProxySubMesh(i);
152     }
153     void setNode2Node(const SMDS_MeshNode*                 srcNode,
154                       const SMDS_MeshNode*                 proxyNode,
155                       const SMESH_ProxyMesh::SubMesh* subMesh)
156     {
157       SMESH_ProxyMesh::setNode2Node( srcNode,proxyNode,subMesh);
158     }
159   };
160   //--------------------------------------------------------------------------------
161   /*!
162    * \brief Listener of events of 3D sub-meshes computed with viscous layers.
163    * It is used to clear an inferior dim sub-meshes modified by viscous layers
164    */
165   class _ShrinkShapeListener : SMESH_subMeshEventListener
166   {
167     _ShrinkShapeListener()
168       : SMESH_subMeshEventListener(/*isDeletable=*/false,
169                                    "StdMeshers_ViscousLayers::_ShrinkShapeListener") {}
170   public:
171     static SMESH_subMeshEventListener* Get() { static _ShrinkShapeListener l; return &l; }
172     virtual void ProcessEvent(const int                       event,
173                               const int                       eventType,
174                               SMESH_subMesh*                  solidSM,
175                               SMESH_subMeshEventListenerData* data,
176                               const SMESH_Hypothesis*         hyp)
177     {
178       if ( SMESH_subMesh::COMPUTE_EVENT == eventType && solidSM->IsEmpty() && data )
179       {
180         SMESH_subMeshEventListener::ProcessEvent(event,eventType,solidSM,data,hyp);
181       }
182     }
183   };
184   //--------------------------------------------------------------------------------
185   /*!
186    * \brief Listener of events of 3D sub-meshes computed with viscous layers.
187    * It is used to store data computed by _ViscousBuilder for a sub-mesh and to
188    * delete the data as soon as it has been used
189    */
190   class _ViscousListener : SMESH_subMeshEventListener
191   {
192     _ViscousListener():
193       SMESH_subMeshEventListener(/*isDeletable=*/false,
194                                  "StdMeshers_ViscousLayers::_ViscousListener") {}
195     static SMESH_subMeshEventListener* Get() { static _ViscousListener l; return &l; }
196   public:
197     virtual void ProcessEvent(const int                       event,
198                               const int                       eventType,
199                               SMESH_subMesh*                  subMesh,
200                               SMESH_subMeshEventListenerData* data,
201                               const SMESH_Hypothesis*         hyp)
202     {
203       if (( SMESH_subMesh::COMPUTE_EVENT       == eventType ) &&
204           ( SMESH_subMesh::CHECK_COMPUTE_STATE != event &&
205             SMESH_subMesh::SUBMESH_COMPUTED    != event ))
206       {
207         // delete SMESH_ProxyMesh containing temporary faces
208         subMesh->DeleteEventListener( this );
209       }
210     }
211     // Finds or creates proxy mesh of the solid
212     static _MeshOfSolid* GetSolidMesh(SMESH_Mesh*         mesh,
213                                       const TopoDS_Shape& solid,
214                                       bool                toCreate=false)
215     {
216       if ( !mesh ) return 0;
217       SMESH_subMesh* sm = mesh->GetSubMesh(solid);
218       _MeshOfSolid* data = (_MeshOfSolid*) sm->GetEventListenerData( Get() );
219       if ( !data && toCreate )
220       {
221         data = new _MeshOfSolid(mesh);
222         data->mySubMeshes.push_back( sm ); // to find SOLID by _MeshOfSolid
223         sm->SetEventListener( Get(), data, sm );
224       }
225       return data;
226     }
227     // Removes proxy mesh of the solid
228     static void RemoveSolidMesh(SMESH_Mesh* mesh, const TopoDS_Shape& solid)
229     {
230       mesh->GetSubMesh(solid)->DeleteEventListener( _ViscousListener::Get() );
231     }
232   };
233   
234   //================================================================================
235   /*!
236    * \brief sets a sub-mesh event listener to clear sub-meshes of sub-shapes of
237    * the main shape when sub-mesh of the main shape is cleared,
238    * for example to clear sub-meshes of FACEs when sub-mesh of a SOLID
239    * is cleared
240    */
241   //================================================================================
242
243   void ToClearSubWithMain( SMESH_subMesh* sub, const TopoDS_Shape& main)
244   {
245     SMESH_subMesh* mainSM = sub->GetFather()->GetSubMesh( main );
246     SMESH_subMeshEventListenerData* data =
247       mainSM->GetEventListenerData( _ShrinkShapeListener::Get());
248     if ( data )
249     {
250       if ( find( data->mySubMeshes.begin(), data->mySubMeshes.end(), sub ) ==
251            data->mySubMeshes.end())
252         data->mySubMeshes.push_back( sub );
253     }
254     else
255     {
256       data = SMESH_subMeshEventListenerData::MakeData( /*dependent=*/sub );
257       sub->SetEventListener( _ShrinkShapeListener::Get(), data, /*whereToListenTo=*/mainSM );
258     }
259   }
260   struct _SolidData;
261   //--------------------------------------------------------------------------------
262   /*!
263    * \brief Simplex (triangle or tetrahedron) based on 1 (tria) or 2 (tet) nodes of
264    * _LayerEdge and 2 nodes of the mesh surface beening smoothed.
265    * The class is used to check validity of face or volumes around a smoothed node;
266    * it stores only 2 nodes as the other nodes are stored by _LayerEdge.
267    */
268   struct _Simplex
269   {
270     const SMDS_MeshNode *_nPrev, *_nNext; // nodes on a smoothed mesh surface
271     const SMDS_MeshNode *_nOpp; // in 2D case, a node opposite to a smoothed node in QUAD
272     _Simplex(const SMDS_MeshNode* nPrev=0,
273              const SMDS_MeshNode* nNext=0,
274              const SMDS_MeshNode* nOpp=0)
275       : _nPrev(nPrev), _nNext(nNext), _nOpp(nOpp) {}
276     bool IsForward(const gp_XYZ* pntSrc, const gp_XYZ* pntTgt, double& vol) const
277     {
278       const double M[3][3] =
279         {{ _nNext->X() - pntSrc->X(), _nNext->Y() - pntSrc->Y(), _nNext->Z() - pntSrc->Z() },
280          { pntTgt->X() - pntSrc->X(), pntTgt->Y() - pntSrc->Y(), pntTgt->Z() - pntSrc->Z() },
281          { _nPrev->X() - pntSrc->X(), _nPrev->Y() - pntSrc->Y(), _nPrev->Z() - pntSrc->Z() }};
282       vol = ( + M[0][0] * M[1][1] * M[2][2]
283               + M[0][1] * M[1][2] * M[2][0]
284               + M[0][2] * M[1][0] * M[2][1]
285               - M[0][0] * M[1][2] * M[2][1]
286               - M[0][1] * M[1][0] * M[2][2]
287               - M[0][2] * M[1][1] * M[2][0]);
288       return vol > 1e-100;
289     }
290     bool IsForward(const SMDS_MeshNode* nSrc, const gp_XYZ& pTgt, double& vol) const
291     {
292       SMESH_TNodeXYZ pSrc( nSrc );
293       return IsForward( &pSrc, &pTgt, vol );
294     }
295     bool IsForward(const gp_XY&         tgtUV,
296                    const SMDS_MeshNode* smoothedNode,
297                    const TopoDS_Face&   face,
298                    SMESH_MesherHelper&  helper,
299                    const double         refSign) const
300     {
301       gp_XY prevUV = helper.GetNodeUV( face, _nPrev, smoothedNode );
302       gp_XY nextUV = helper.GetNodeUV( face, _nNext, smoothedNode );
303       gp_Vec2d v1( tgtUV, prevUV ), v2( tgtUV, nextUV );
304       double d = v1 ^ v2;
305       return d*refSign > 1e-100;
306     }
307     bool IsMinAngleOK( const gp_XYZ& pTgt, double& minAngle ) const
308     {
309       SMESH_TNodeXYZ pPrev( _nPrev ), pNext( _nNext );
310       if ( !_nOpp ) // triangle
311       {
312         gp_Vec tp( pPrev - pTgt ), pn( pNext - pPrev ), nt( pTgt - pNext );
313         double tp2 = tp.SquareMagnitude();
314         double pn2 = pn.SquareMagnitude();
315         double nt2 = nt.SquareMagnitude();
316
317         if ( tp2 < pn2 && tp2 < nt2 )
318           minAngle = ( nt * -pn ) * ( nt * -pn ) / nt2 / pn2;
319         else if ( pn2 < nt2 )
320           minAngle = ( tp * -nt ) * ( tp * -nt ) / tp2 / nt2;
321         else
322           minAngle = ( pn * -tp ) * ( pn * -tp ) / pn2 / tp2;
323
324         static double theMaxCos2 = ( Cos( theMinSmoothTriaAngle * M_PI / 180. ) *
325                                      Cos( theMinSmoothTriaAngle * M_PI / 180. ));
326         return minAngle < theMaxCos2;
327       }
328       else // quadrangle
329       {
330         SMESH_TNodeXYZ pOpp( _nOpp );
331         gp_Vec tp( pPrev - pTgt ), po( pOpp - pPrev ), on( pNext - pOpp), nt( pTgt - pNext );
332         double tp2 = tp.SquareMagnitude();
333         double po2 = po.SquareMagnitude();
334         double on2 = on.SquareMagnitude();
335         double nt2 = nt.SquareMagnitude();
336         minAngle = Max( Max((( tp * -nt ) * ( tp * -nt ) / tp2 / nt2 ),
337                             (( po * -tp ) * ( po * -tp ) / po2 / tp2 )),
338                         Max((( on * -po ) * ( on * -po ) / on2 / po2 ),
339                             (( nt * -on ) * ( nt * -on ) / nt2 / on2 )));
340
341         static double theMaxCos2 = ( Cos( theMinSmoothQuadAngle * M_PI / 180. ) *
342                                      Cos( theMinSmoothQuadAngle * M_PI / 180. ));
343         return minAngle < theMaxCos2;
344       }
345     }
346     bool IsNeighbour(const _Simplex& other) const
347     {
348       return _nPrev == other._nNext || _nNext == other._nPrev;
349     }
350     bool Includes( const SMDS_MeshNode* node ) const { return _nPrev == node || _nNext == node; }
351     static void GetSimplices( const SMDS_MeshNode* node,
352                               vector<_Simplex>&   simplices,
353                               const set<TGeomID>& ingnoreShapes,
354                               const _SolidData*   dataToCheckOri = 0,
355                               const bool          toSort = false);
356     static void SortSimplices(vector<_Simplex>& simplices);
357   };
358   //--------------------------------------------------------------------------------
359   /*!
360    * Structure used to take into account surface curvature while smoothing
361    */
362   struct _Curvature
363   {
364     double   _r; // radius
365     double   _k; // factor to correct node smoothed position
366     double   _h2lenRatio; // avgNormProj / (2*avgDist)
367     gp_Pnt2d _uv; // UV used in putOnOffsetSurface()
368   public:
369     static _Curvature* New( double avgNormProj, double avgDist )
370     {
371       _Curvature* c = 0;
372       if ( fabs( avgNormProj / avgDist ) > 1./200 )
373       {
374         c = new _Curvature;
375         c->_r = avgDist * avgDist / avgNormProj;
376         c->_k = avgDist * avgDist / c->_r / c->_r;
377         //c->_k = avgNormProj / c->_r;
378         c->_k *= ( c->_r < 0 ? 1/1.1 : 1.1 ); // not to be too restrictive
379         c->_h2lenRatio = avgNormProj / ( avgDist + avgDist );
380
381         c->_uv.SetCoord( 0., 0. );
382       }
383       return c;
384     }
385     double lenDelta(double len) const { return _k * ( _r + len ); }
386     double lenDeltaByDist(double dist) const { return dist * _h2lenRatio; }
387   };
388   //--------------------------------------------------------------------------------
389
390   struct _2NearEdges;
391   struct _LayerEdge;
392   struct _EdgesOnShape;
393   struct _Smoother1D;
394   typedef map< const SMDS_MeshNode*, _LayerEdge*, TIDCompare > TNode2Edge;
395
396   //--------------------------------------------------------------------------------
397   /*!
398    * \brief Edge normal to surface, connecting a node on solid surface (_nodes[0])
399    * and a node of the most internal layer (_nodes.back())
400    */
401   struct _LayerEdge
402   {
403     typedef gp_XYZ (_LayerEdge::*PSmooFun)();
404
405     vector< const SMDS_MeshNode*> _nodes;
406
407     gp_XYZ              _normal;    // to boundary of solid
408     vector<gp_XYZ>      _pos;       // points computed during inflation
409     double              _len;       // length achived with the last inflation step
410     double              _maxLen;    // maximal possible length
411     double              _cosin;     // of angle (_normal ^ surface)
412     double              _minAngle;  // of _simplices
413     double              _lenFactor; // to compute _len taking _cosin into account
414     int                 _flags;
415
416     // simplices connected to the source node (_nodes[0]);
417     // used for smoothing and quality check of _LayerEdge's based on the FACE
418     vector<_Simplex>    _simplices;
419     vector<_LayerEdge*> _neibors; // all surrounding _LayerEdge's
420     PSmooFun            _smooFunction; // smoothing function
421     _Curvature*         _curvature;
422     // data for smoothing of _LayerEdge's based on the EDGE
423     _2NearEdges*        _2neibors;
424
425     enum EFlags { TO_SMOOTH       = 1,
426                   MOVED           = 2,   // set by _neibors[i]->SetNewLength()
427                   SMOOTHED        = 4,   // set by this->Smooth()
428                   DIFFICULT       = 8,   // near concave VERTEX
429                   ON_CONCAVE_FACE = 16,
430                   BLOCKED         = 32,  // not to inflate any more
431                   INTERSECTED     = 64,  // close intersection with a face found
432                   NORMAL_UPDATED  = 128,
433                   MARKED          = 256, // local usage
434                   MULTI_NORMAL    = 512, // a normal is invisible by some of surrounding faces
435                   NEAR_BOUNDARY   = 1024,// is near FACE boundary forcing smooth
436                   SMOOTHED_C1     = 2048,// is on _eosC1
437                   DISTORTED       = 4096,// was bad before smoothing
438                   RISKY_SWOL      = 8192 // SWOL is parallel to a source FACE
439     };
440     bool Is   ( EFlags f ) const { return _flags & f; }
441     void Set  ( EFlags f ) { _flags |= f; }
442     void Unset( EFlags f ) { _flags &= ~f; }
443
444     void SetNewLength( double len, _EdgesOnShape& eos, SMESH_MesherHelper& helper );
445     bool SetNewLength2d( Handle(Geom_Surface)& surface,
446                          const TopoDS_Face&    F,
447                          _EdgesOnShape&        eos,
448                          SMESH_MesherHelper&   helper );
449     void SetDataByNeighbors( const SMDS_MeshNode* n1,
450                              const SMDS_MeshNode* n2,
451                              const _EdgesOnShape& eos,
452                              SMESH_MesherHelper&  helper);
453     void Block( _SolidData& data );
454     void InvalidateStep( size_t curStep, const _EdgesOnShape& eos, bool restoreLength=false );
455     void ChooseSmooFunction(const set< TGeomID >& concaveVertices,
456                             const TNode2Edge&     n2eMap);
457     void SmoothPos( const vector< double >& segLen, const double tol );
458     int  Smooth(const int step, const bool isConcaveFace, bool findBest);
459     int  Smooth(const int step, bool findBest, vector< _LayerEdge* >& toSmooth );
460     int  CheckNeiborsOnBoundary(vector< _LayerEdge* >* badNeibors = 0 );
461     void SmoothWoCheck();
462     bool SmoothOnEdge(Handle(ShapeAnalysis_Surface)& surface,
463                       const TopoDS_Face&             F,
464                       SMESH_MesherHelper&            helper);
465     void MoveNearConcaVer( const _EdgesOnShape*    eov,
466                            const _EdgesOnShape*    eos,
467                            const int               step,
468                            vector< _LayerEdge* > & badSmooEdges);
469     bool FindIntersection( SMESH_ElementSearcher&   searcher,
470                            double &                 distance,
471                            const double&            epsilon,
472                            _EdgesOnShape&           eos,
473                            const SMDS_MeshElement** face = 0);
474     bool SegTriaInter( const gp_Ax1&        lastSegment,
475                        const gp_XYZ&        p0,
476                        const gp_XYZ&        p1,
477                        const gp_XYZ&        p2,
478                        double&              dist,
479                        const double&        epsilon) const;
480     bool SegTriaInter( const gp_Ax1&        lastSegment,
481                        const SMDS_MeshNode* n0,
482                        const SMDS_MeshNode* n1,
483                        const SMDS_MeshNode* n2,
484                        double&              dist,
485                        const double&        epsilon) const
486     { return SegTriaInter( lastSegment,
487                            SMESH_TNodeXYZ( n0 ), SMESH_TNodeXYZ( n1 ), SMESH_TNodeXYZ( n2 ),
488                            dist, epsilon );
489     }
490     const gp_XYZ& PrevPos() const { return _pos[ _pos.size() - 2 ]; }
491     const gp_XYZ& PrevCheckPos() const { return _pos[ Is( NORMAL_UPDATED ) ? _pos.size()-2 : 0 ]; }
492     gp_Ax1 LastSegment(double& segLen, _EdgesOnShape& eos) const;
493     gp_XY  LastUV( const TopoDS_Face& F, _EdgesOnShape& eos ) const;
494     bool   IsOnEdge() const { return _2neibors; }
495     gp_XYZ Copy( _LayerEdge& other, _EdgesOnShape& eos, SMESH_MesherHelper& helper );
496     void   SetCosin( double cosin );
497     void   SetNormal( const gp_XYZ& n ) { _normal = n; }
498     int    NbSteps() const { return _pos.size() - 1; } // nb inlation steps
499     bool   IsNeiborOnEdge( const _LayerEdge* edge ) const;
500     void   SetSmooLen( double len ) { // set _len at which smoothing is needed
501       _cosin = len; // as for _LayerEdge's on FACE _cosin is not used
502     }
503     double GetSmooLen() { return _cosin; } // for _LayerEdge's on FACE _cosin is not used
504
505     gp_XYZ smoothLaplacian();
506     gp_XYZ smoothAngular();
507     gp_XYZ smoothLengthWeighted();
508     gp_XYZ smoothCentroidal();
509     gp_XYZ smoothNefPolygon();
510
511     enum { FUN_LAPLACIAN, FUN_LENWEIGHTED, FUN_CENTROIDAL, FUN_NEFPOLY, FUN_ANGULAR, FUN_NB };
512     static const int theNbSmooFuns = FUN_NB;
513     static PSmooFun _funs[theNbSmooFuns];
514     static const char* _funNames[theNbSmooFuns+1];
515     int smooFunID( PSmooFun fun=0) const;
516   };
517   _LayerEdge::PSmooFun _LayerEdge::_funs[theNbSmooFuns] = { &_LayerEdge::smoothLaplacian,
518                                                             &_LayerEdge::smoothLengthWeighted,
519                                                             &_LayerEdge::smoothCentroidal,
520                                                             &_LayerEdge::smoothNefPolygon,
521                                                             &_LayerEdge::smoothAngular };
522   const char* _LayerEdge::_funNames[theNbSmooFuns+1] = { "Laplacian",
523                                                          "LengthWeighted",
524                                                          "Centroidal",
525                                                          "NefPolygon",
526                                                          "Angular",
527                                                          "None"};
528   struct _LayerEdgeCmp
529   {
530     bool operator () (const _LayerEdge* e1, const _LayerEdge* e2) const
531     {
532       const bool cmpNodes = ( e1 && e2 && e1->_nodes.size() && e2->_nodes.size() );
533       return cmpNodes ? ( e1->_nodes[0]->GetID() < e2->_nodes[0]->GetID()) : ( e1 < e2 );
534     }
535   };
536   //--------------------------------------------------------------------------------
537   /*!
538    * A 2D half plane used by _LayerEdge::smoothNefPolygon()
539    */
540   struct _halfPlane
541   {
542     gp_XY _pos, _dir, _inNorm;
543     bool IsOut( const gp_XY p, const double tol ) const
544     {
545       return _inNorm * ( p - _pos ) < -tol;
546     }
547     bool FindIntersection( const _halfPlane& hp, gp_XY & intPnt )
548     {
549       //const double eps = 1e-10;
550       double D = _dir.Crossed( hp._dir );
551       if ( fabs(D) < std::numeric_limits<double>::min())
552         return false;
553       gp_XY vec21 = _pos - hp._pos; 
554       double u = hp._dir.Crossed( vec21 ) / D; 
555       intPnt = _pos + _dir * u;
556       return true;
557     }
558   };
559   //--------------------------------------------------------------------------------
560   /*!
561    * Structure used to smooth a _LayerEdge based on an EDGE.
562    */
563   struct _2NearEdges
564   {
565     double               _wgt  [2]; // weights of _nodes
566     _LayerEdge*          _edges[2];
567
568      // normal to plane passing through _LayerEdge._normal and tangent of EDGE
569     gp_XYZ*              _plnNorm;
570
571     _2NearEdges() { _edges[0]=_edges[1]=0; _plnNorm = 0; }
572     const SMDS_MeshNode* tgtNode(bool is2nd) {
573       return _edges[is2nd] ? _edges[is2nd]->_nodes.back() : 0;
574     }
575     const SMDS_MeshNode* srcNode(bool is2nd) {
576       return _edges[is2nd] ? _edges[is2nd]->_nodes[0] : 0;
577     }
578     void reverse() {
579       std::swap( _wgt  [0], _wgt  [1] );
580       std::swap( _edges[0], _edges[1] );
581     }
582     void set( _LayerEdge* e1, _LayerEdge* e2, double w1, double w2 ) {
583       _edges[0] = e1; _edges[1] = e2; _wgt[0] = w1; _wgt[1] = w2;
584     }
585     bool include( const _LayerEdge* e ) {
586       return ( _edges[0] == e || _edges[1] == e );
587     }
588   };
589
590
591   //--------------------------------------------------------------------------------
592   /*!
593    * \brief Layers parameters got by averaging several hypotheses
594    */
595   struct AverageHyp
596   {
597     AverageHyp( const StdMeshers_ViscousLayers* hyp = 0 )
598       :_nbLayers(0), _nbHyps(0), _method(0), _thickness(0), _stretchFactor(0)
599     {
600       Add( hyp );
601     }
602     void Add( const StdMeshers_ViscousLayers* hyp )
603     {
604       if ( hyp )
605       {
606         _nbHyps++;
607         _nbLayers       = hyp->GetNumberLayers();
608         //_thickness     += hyp->GetTotalThickness();
609         _thickness      = Max( _thickness, hyp->GetTotalThickness() );
610         _stretchFactor += hyp->GetStretchFactor();
611         _method         = hyp->GetMethod();
612       }
613     }
614     double GetTotalThickness() const { return _thickness; /*_nbHyps ? _thickness / _nbHyps : 0;*/ }
615     double GetStretchFactor()  const { return _nbHyps ? _stretchFactor / _nbHyps : 0; }
616     int    GetNumberLayers()   const { return _nbLayers; }
617     int    GetMethod()         const { return _method; }
618
619     bool   UseSurfaceNormal()  const
620     { return _method == StdMeshers_ViscousLayers::SURF_OFFSET_SMOOTH; }
621     bool   ToSmooth()          const
622     { return _method == StdMeshers_ViscousLayers::SURF_OFFSET_SMOOTH; }
623     bool   IsOffsetMethod()    const
624     { return _method == StdMeshers_ViscousLayers::FACE_OFFSET; }
625
626   private:
627     int     _nbLayers, _nbHyps, _method;
628     double  _thickness, _stretchFactor;
629   };
630
631   //--------------------------------------------------------------------------------
632   /*!
633    * \brief _LayerEdge's on a shape and other shape data
634    */
635   struct _EdgesOnShape
636   {
637     vector< _LayerEdge* > _edges;
638
639     TopoDS_Shape          _shape;
640     TGeomID               _shapeID;
641     SMESH_subMesh *       _subMesh;
642     // face or edge w/o layer along or near which _edges are inflated
643     TopoDS_Shape          _sWOL;
644     bool                  _isRegularSWOL; // w/o singularities
645     // averaged StdMeshers_ViscousLayers parameters
646     AverageHyp            _hyp;
647     bool                  _toSmooth;
648     _Smoother1D*          _edgeSmoother;
649     vector< _EdgesOnShape* > _eosConcaVer; // edges at concave VERTEXes of a FACE
650     vector< _EdgesOnShape* > _eosC1; // to smooth together several C1 continues shapes
651
652     vector< gp_XYZ >         _faceNormals; // if _shape is FACE
653     vector< _EdgesOnShape* > _faceEOS; // to get _faceNormals of adjacent FACEs
654
655     Handle(ShapeAnalysis_Surface) _offsetSurf;
656     _LayerEdge*                   _edgeForOffset;
657
658     _SolidData*            _data; // parent SOLID
659
660     TopAbs_ShapeEnum ShapeType() const
661     { return _shape.IsNull() ? TopAbs_SHAPE : _shape.ShapeType(); }
662     TopAbs_ShapeEnum SWOLType() const
663     { return _sWOL.IsNull() ? TopAbs_SHAPE : _sWOL.ShapeType(); }
664     bool             HasC1( const _EdgesOnShape* other ) const
665     { return std::find( _eosC1.begin(), _eosC1.end(), other ) != _eosC1.end(); }
666     bool             GetNormal( const SMDS_MeshElement* face, gp_Vec& norm );
667     _SolidData&      GetData() const { return *_data; }
668
669     _EdgesOnShape(): _shapeID(-1), _subMesh(0), _toSmooth(false), _edgeSmoother(0) {}
670   };
671
672   //--------------------------------------------------------------------------------
673   /*!
674    * \brief Convex FACE whose radius of curvature is less than the thickness of 
675    *        layers. It is used to detect distortion of prisms based on a convex
676    *        FACE and to update normals to enable further increasing the thickness
677    */
678   struct _ConvexFace
679   {
680     TopoDS_Face                     _face;
681
682     // edges whose _simplices are used to detect prism distortion
683     vector< _LayerEdge* >           _simplexTestEdges;
684
685     // map a sub-shape to _SolidData::_edgesOnShape
686     map< TGeomID, _EdgesOnShape* >  _subIdToEOS;
687
688     bool                            _normalsFixed;
689
690     bool GetCenterOfCurvature( _LayerEdge*         ledge,
691                                BRepLProp_SLProps&  surfProp,
692                                SMESH_MesherHelper& helper,
693                                gp_Pnt &            center ) const;
694     bool CheckPrisms() const;
695   };
696
697   //--------------------------------------------------------------------------------
698   /*!
699    * \brief Structure holding _LayerEdge's based on EDGEs that will collide
700    *        at inflation up to the full thickness. A detected collision
701    *        is fixed in updateNormals()
702    */
703   struct _CollisionEdges
704   {
705     _LayerEdge*           _edge;
706     vector< _LayerEdge* > _intEdges; // each pair forms an intersected quadrangle
707     const SMDS_MeshNode* nSrc(int i) const { return _intEdges[i]->_nodes[0]; }
708     const SMDS_MeshNode* nTgt(int i) const { return _intEdges[i]->_nodes.back(); }
709   };
710
711   //--------------------------------------------------------------------------------
712   /*!
713    * \brief Data of a SOLID
714    */
715   struct _SolidData
716   {
717     typedef const StdMeshers_ViscousLayers* THyp;
718     TopoDS_Shape                    _solid;
719     TGeomID                         _index; // SOLID id
720     _MeshOfSolid*                   _proxyMesh;
721     list< THyp >                    _hyps;
722     list< TopoDS_Shape >            _hypShapes;
723     map< TGeomID, THyp >            _face2hyp; // filled if _hyps.size() > 1
724     set< TGeomID >                  _reversedFaceIds;
725     set< TGeomID >                  _ignoreFaceIds; // WOL FACEs and FACEs of other SOLIDs
726
727     double                          _stepSize, _stepSizeCoeff, _geomSize;
728     const SMDS_MeshNode*            _stepSizeNodes[2];
729
730     TNode2Edge                      _n2eMap; // nodes and _LayerEdge's based on them
731
732     // map to find _n2eMap of another _SolidData by a shrink shape shared by two _SolidData's
733     map< TGeomID, TNode2Edge* >     _s2neMap;
734     // _LayerEdge's with underlying shapes
735     vector< _EdgesOnShape >         _edgesOnShape;
736
737     // key:   an id of shape (EDGE or VERTEX) shared by a FACE with
738     //        layers and a FACE w/o layers
739     // value: the shape (FACE or EDGE) to shrink mesh on.
740     //       _LayerEdge's basing on nodes on key shape are inflated along the value shape
741     map< TGeomID, TopoDS_Shape >     _shrinkShape2Shape;
742
743     // Convex FACEs whose radius of curvature is less than the thickness of layers
744     map< TGeomID, _ConvexFace >      _convexFaces;
745
746     // shapes (EDGEs and VERTEXes) srink from which is forbidden due to collisions with
747     // the adjacent SOLID
748     set< TGeomID >                   _noShrinkShapes;
749
750     int                              _nbShapesToSmooth;
751
752     //map< TGeomID,Handle(Geom_Curve)> _edge2curve;
753
754     vector< _CollisionEdges >        _collisionEdges;
755     set< TGeomID >                   _concaveFaces;
756
757     double                           _maxThickness; // of all _hyps
758     double                           _minThickness; // of all _hyps
759
760     double                           _epsilon; // precision for SegTriaInter()
761
762     SMESH_MesherHelper*              _helper;
763
764     _SolidData(const TopoDS_Shape& s=TopoDS_Shape(),
765                _MeshOfSolid*       m=0)
766       :_solid(s), _proxyMesh(m), _helper(0) {}
767     ~_SolidData();
768
769     void SortOnEdge( const TopoDS_Edge& E, vector< _LayerEdge* >& edges);
770     void Sort2NeiborsOnEdge( vector< _LayerEdge* >& edges );
771
772     _ConvexFace* GetConvexFace( const TGeomID faceID ) {
773       map< TGeomID, _ConvexFace >::iterator id2face = _convexFaces.find( faceID );
774       return id2face == _convexFaces.end() ? 0 : & id2face->second;
775     }
776     _EdgesOnShape* GetShapeEdges(const TGeomID       shapeID );
777     _EdgesOnShape* GetShapeEdges(const TopoDS_Shape& shape );
778     _EdgesOnShape* GetShapeEdges(const _LayerEdge*   edge )
779     { return GetShapeEdges( edge->_nodes[0]->getshapeId() ); }
780
781     SMESH_MesherHelper& GetHelper() const { return *_helper; }
782
783     void UnmarkEdges() {
784       for ( size_t i = 0; i < _edgesOnShape.size(); ++i )
785         for ( size_t j = 0; j < _edgesOnShape[i]._edges.size(); ++j )
786           _edgesOnShape[i]._edges[j]->Unset( _LayerEdge::MARKED );
787     }
788     void AddShapesToSmooth( const set< _EdgesOnShape* >& shape,
789                             const set< _EdgesOnShape* >* edgesNoAnaSmooth=0 );
790
791     void PrepareEdgesToSmoothOnFace( _EdgesOnShape* eof, bool substituteSrcNodes );
792   };
793   //--------------------------------------------------------------------------------
794   /*!
795    * \brief Offset plane used in getNormalByOffset()
796    */
797   struct _OffsetPlane
798   {
799     gp_Pln _plane;
800     int    _faceIndex;
801     int    _faceIndexNext[2];
802     gp_Lin _lines[2]; // line of intersection with neighbor _OffsetPlane's
803     bool   _isLineOK[2];
804     _OffsetPlane() {
805       _isLineOK[0] = _isLineOK[1] = false; _faceIndexNext[0] = _faceIndexNext[1] = -1;
806     }
807     void   ComputeIntersectionLine( _OffsetPlane& pln );
808     gp_XYZ GetCommonPoint(bool& isFound) const;
809     int    NbLines() const { return _isLineOK[0] + _isLineOK[1]; }
810   };
811   //--------------------------------------------------------------------------------
812   /*!
813    * \brief Container of centers of curvature at nodes on an EDGE bounding _ConvexFace
814    */
815   struct _CentralCurveOnEdge
816   {
817     bool                  _isDegenerated;
818     vector< gp_Pnt >      _curvaCenters;
819     vector< _LayerEdge* > _ledges;
820     vector< gp_XYZ >      _normals; // new normal for each of _ledges
821     vector< double >      _segLength2;
822
823     TopoDS_Edge           _edge;
824     TopoDS_Face           _adjFace;
825     bool                  _adjFaceToSmooth;
826
827     void Append( const gp_Pnt& center, _LayerEdge* ledge )
828     {
829       if ( _curvaCenters.size() > 0 )
830         _segLength2.push_back( center.SquareDistance( _curvaCenters.back() ));
831       _curvaCenters.push_back( center );
832       _ledges.push_back( ledge );
833       _normals.push_back( ledge->_normal );
834     }
835     bool FindNewNormal( const gp_Pnt& center, gp_XYZ& newNormal );
836     void SetShapes( const TopoDS_Edge&  edge,
837                     const _ConvexFace&  convFace,
838                     _SolidData&         data,
839                     SMESH_MesherHelper& helper);
840   };
841   //--------------------------------------------------------------------------------
842   /*!
843    * \brief Data of node on a shrinked FACE
844    */
845   struct _SmoothNode
846   {
847     const SMDS_MeshNode*         _node;
848     vector<_Simplex>             _simplices; // for quality check
849
850     enum SmoothType { LAPLACIAN, CENTROIDAL, ANGULAR, TFI };
851
852     bool Smooth(int&                  badNb,
853                 Handle(Geom_Surface)& surface,
854                 SMESH_MesherHelper&   helper,
855                 const double          refSign,
856                 SmoothType            how,
857                 bool                  set3D);
858
859     gp_XY computeAngularPos(vector<gp_XY>& uv,
860                             const gp_XY&   uvToFix,
861                             const double   refSign );
862   };
863   //--------------------------------------------------------------------------------
864   /*!
865    * \brief Builder of viscous layers
866    */
867   class _ViscousBuilder
868   {
869   public:
870     _ViscousBuilder();
871     // does it's job
872     SMESH_ComputeErrorPtr Compute(SMESH_Mesh&         mesh,
873                                   const TopoDS_Shape& shape);
874     // check validity of hypotheses
875     SMESH_ComputeErrorPtr CheckHypotheses( SMESH_Mesh&         mesh,
876                                            const TopoDS_Shape& shape );
877
878     // restore event listeners used to clear an inferior dim sub-mesh modified by viscous layers
879     void RestoreListeners();
880
881     // computes SMESH_ProxyMesh::SubMesh::_n2n;
882     bool MakeN2NMap( _MeshOfSolid* pm );
883
884   private:
885
886     bool findSolidsWithLayers();
887     bool findFacesWithLayers(const bool onlyWith=false);
888     void getIgnoreFaces(const TopoDS_Shape&             solid,
889                         const StdMeshers_ViscousLayers* hyp,
890                         const TopoDS_Shape&             hypShape,
891                         set<TGeomID>&                   ignoreFaces);
892     bool makeLayer(_SolidData& data);
893     void setShapeData( _EdgesOnShape& eos, SMESH_subMesh* sm, _SolidData& data );
894     bool setEdgeData( _LayerEdge& edge, _EdgesOnShape& eos,
895                       SMESH_MesherHelper& helper, _SolidData& data);
896     gp_XYZ getFaceNormal(const SMDS_MeshNode* n,
897                          const TopoDS_Face&   face,
898                          SMESH_MesherHelper&  helper,
899                          bool&                isOK,
900                          bool                 shiftInside=false);
901     bool getFaceNormalAtSingularity(const gp_XY&        uv,
902                                     const TopoDS_Face&  face,
903                                     SMESH_MesherHelper& helper,
904                                     gp_Dir&             normal );
905     gp_XYZ getWeigthedNormal( const _LayerEdge*                edge );
906     gp_XYZ getNormalByOffset( _LayerEdge*                      edge,
907                               std::pair< TopoDS_Face, gp_XYZ > fId2Normal[],
908                               int                              nbFaces );
909     bool findNeiborsOnEdge(const _LayerEdge*     edge,
910                            const SMDS_MeshNode*& n1,
911                            const SMDS_MeshNode*& n2,
912                            _EdgesOnShape&        eos,
913                            _SolidData&           data);
914     void findSimplexTestEdges( _SolidData&                    data,
915                                vector< vector<_LayerEdge*> >& edgesByGeom);
916     void computeGeomSize( _SolidData& data );
917     bool findShapesToSmooth( _SolidData& data);
918     void limitStepSizeByCurvature( _SolidData&  data );
919     void limitStepSize( _SolidData&             data,
920                         const SMDS_MeshElement* face,
921                         const _LayerEdge*       maxCosinEdge );
922     void limitStepSize( _SolidData& data, const double minSize);
923     bool inflate(_SolidData& data);
924     bool smoothAndCheck(_SolidData& data, const int nbSteps, double & distToIntersection);
925     int  invalidateBadSmooth( _SolidData&               data,
926                               SMESH_MesherHelper&       helper,
927                               vector< _LayerEdge* >&    badSmooEdges,
928                               vector< _EdgesOnShape* >& eosC1,
929                               const int                 infStep );
930     void makeOffsetSurface( _EdgesOnShape& eos, SMESH_MesherHelper& );
931     void putOnOffsetSurface( _EdgesOnShape& eos, int infStep, int smooStep=0, bool moveAll=false );
932     void findCollisionEdges( _SolidData& data, SMESH_MesherHelper& helper );
933     bool updateNormals( _SolidData& data, SMESH_MesherHelper& helper, int stepNb, double stepSize );
934     bool updateNormalsOfConvexFaces( _SolidData&         data,
935                                      SMESH_MesherHelper& helper,
936                                      int                 stepNb );
937     void updateNormalsOfC1Vertices( _SolidData& data );
938     bool updateNormalsOfSmoothed( _SolidData&         data,
939                                   SMESH_MesherHelper& helper,
940                                   const int           nbSteps,
941                                   const double        stepSize );
942     bool isNewNormalOk( _SolidData&   data,
943                         _LayerEdge&   edge,
944                         const gp_XYZ& newNormal);
945     bool refine(_SolidData& data);
946     bool shrink();
947     bool prepareEdgeToShrink( _LayerEdge& edge, _EdgesOnShape& eos,
948                               SMESH_MesherHelper& helper,
949                               const SMESHDS_SubMesh* faceSubMesh );
950     void restoreNoShrink( _LayerEdge& edge ) const;
951     void fixBadFaces(const TopoDS_Face&          F,
952                      SMESH_MesherHelper&         helper,
953                      const bool                  is2D,
954                      const int                   step,
955                      set<const SMDS_MeshNode*> * involvedNodes=NULL);
956     bool addBoundaryElements();
957
958     bool error( const string& text, int solidID=-1 );
959     SMESHDS_Mesh* getMeshDS() const { return _mesh->GetMeshDS(); }
960
961     // debug
962     void makeGroupOfLE();
963
964     SMESH_Mesh*           _mesh;
965     SMESH_ComputeErrorPtr _error;
966
967     vector< _SolidData >  _sdVec;
968     int                   _tmpFaceID;
969   };
970   //--------------------------------------------------------------------------------
971   /*!
972    * \brief Shrinker of nodes on the EDGE
973    */
974   class _Shrinker1D
975   {
976     TopoDS_Edge                   _geomEdge;
977     vector<double>                _initU;
978     vector<double>                _normPar;
979     vector<const SMDS_MeshNode*>  _nodes;
980     const _LayerEdge*             _edges[2];
981     bool                          _done;
982   public:
983     void AddEdge( const _LayerEdge* e, _EdgesOnShape& eos, SMESH_MesherHelper& helper );
984     void Compute(bool set3D, SMESH_MesherHelper& helper);
985     void RestoreParams();
986     void SwapSrcTgtNodes(SMESHDS_Mesh* mesh);
987     const TopoDS_Edge& GeomEdge() const { return _geomEdge; }
988     const SMDS_MeshNode* TgtNode( bool is2nd ) const
989     { return _edges[is2nd] ? _edges[is2nd]->_nodes.back() : 0; }
990     const SMDS_MeshNode* SrcNode( bool is2nd ) const
991     { return _edges[is2nd] ? _edges[is2nd]->_nodes[0] : 0; }
992   };
993   //--------------------------------------------------------------------------------
994   /*!
995    * \brief Smoother of _LayerEdge's on EDGE.
996    */
997   struct _Smoother1D
998   {
999     struct OffPnt // point of the offsetted EDGE
1000     {
1001       gp_XYZ      _xyz;    // coord of a point inflated from EDGE w/o smooth
1002       double      _len;    // length reached at previous inflation step
1003       _2NearEdges _2edges; // 2 neighbor _LayerEdge's
1004       double Distance( const OffPnt& p ) const { return ( _xyz - p._xyz ).Modulus(); }
1005     };
1006     vector< OffPnt >   _offPoints;
1007     vector< double >   _leParams; // normalized param of _eos._edges on EDGE
1008     Handle(Geom_Curve) _anaCurve; // for analytic smooth
1009     _LayerEdge         _leOnV[2]; // _LayerEdge's holding normal to the EDGE at VERTEXes
1010     size_t             _iSeg[2];  // index of segment where extreme tgt node is projected
1011     _EdgesOnShape&     _eos;
1012
1013     static Handle(Geom_Curve) CurveForSmooth( const TopoDS_Edge&  E,
1014                                               _EdgesOnShape&      eos,
1015                                               SMESH_MesherHelper& helper);
1016
1017     _Smoother1D( Handle(Geom_Curve) curveForSmooth,
1018                  _EdgesOnShape&     eos )
1019       : _anaCurve( curveForSmooth ), _eos( eos )
1020     {
1021     }
1022     bool Perform(_SolidData&                    data,
1023                  Handle(ShapeAnalysis_Surface)& surface,
1024                  const TopoDS_Face&             F,
1025                  SMESH_MesherHelper&            helper )
1026     {
1027       if ( _leParams.empty() || ( !isAnalytic() && _offPoints.empty() ))
1028         prepare( data );
1029
1030       if ( isAnalytic() )
1031         return smoothAnalyticEdge( data, surface, F, helper );
1032       else
1033         return smoothComplexEdge ( data, surface, F, helper );
1034     }
1035     void prepare(_SolidData& data );
1036
1037     bool smoothAnalyticEdge( _SolidData&                    data,
1038                              Handle(ShapeAnalysis_Surface)& surface,
1039                              const TopoDS_Face&             F,
1040                              SMESH_MesherHelper&            helper);
1041
1042     bool smoothComplexEdge( _SolidData&                    data,
1043                             Handle(ShapeAnalysis_Surface)& surface,
1044                             const TopoDS_Face&             F,
1045                             SMESH_MesherHelper&            helper);
1046
1047     void setNormalOnV( const bool          is2nd,
1048                        SMESH_MesherHelper& helper);
1049
1050     _LayerEdge* getLEdgeOnV( bool is2nd )
1051     {
1052       return _eos._edges[ is2nd ? _eos._edges.size()-1 : 0 ]->_2neibors->_edges[ is2nd ];
1053     }
1054     bool isAnalytic() const { return !_anaCurve.IsNull(); }
1055   };
1056   //--------------------------------------------------------------------------------
1057   /*!
1058    * \brief Class of temporary mesh face.
1059    * We can't use SMDS_FaceOfNodes since it's impossible to set it's ID which is
1060    * needed because SMESH_ElementSearcher internaly uses set of elements sorted by ID
1061    */
1062   struct _TmpMeshFace : public SMDS_MeshElement
1063   {
1064     vector<const SMDS_MeshNode* > _nn;
1065     _TmpMeshFace( const vector<const SMDS_MeshNode*>& nodes,
1066                   int id, int faceID=-1, int idInFace=-1):
1067       SMDS_MeshElement(id), _nn(nodes) { setShapeId(faceID); setIdInShape(idInFace); }
1068     virtual const SMDS_MeshNode* GetNode(const int ind) const { return _nn[ind]; }
1069     virtual SMDSAbs_ElementType  GetType() const              { return SMDSAbs_Face; }
1070     virtual vtkIdType GetVtkType() const                      { return -1; }
1071     virtual SMDSAbs_EntityType   GetEntityType() const        { return SMDSEntity_Last; }
1072     virtual SMDSAbs_GeometryType GetGeomType() const
1073     { return _nn.size() == 3 ? SMDSGeom_TRIANGLE : SMDSGeom_QUADRANGLE; }
1074     virtual SMDS_ElemIteratorPtr elementsIterator(SMDSAbs_ElementType) const
1075     { return SMDS_ElemIteratorPtr( new SMDS_NodeVectorElemIterator( _nn.begin(), _nn.end()));}
1076   };
1077   //--------------------------------------------------------------------------------
1078   /*!
1079    * \brief Class of temporary mesh face storing _LayerEdge it's based on
1080    */
1081   struct _TmpMeshFaceOnEdge : public _TmpMeshFace
1082   {
1083     _LayerEdge *_le1, *_le2;
1084     _TmpMeshFaceOnEdge( _LayerEdge* le1, _LayerEdge* le2, int ID ):
1085       _TmpMeshFace( vector<const SMDS_MeshNode*>(4), ID ), _le1(le1), _le2(le2)
1086     {
1087       _nn[0]=_le1->_nodes[0];
1088       _nn[1]=_le1->_nodes.back();
1089       _nn[2]=_le2->_nodes.back();
1090       _nn[3]=_le2->_nodes[0];
1091     }
1092     gp_XYZ GetDir() const // return average direction of _LayerEdge's, normal to EDGE
1093     {
1094       SMESH_TNodeXYZ p0s( _nn[0] );
1095       SMESH_TNodeXYZ p0t( _nn[1] );
1096       SMESH_TNodeXYZ p1t( _nn[2] );
1097       SMESH_TNodeXYZ p1s( _nn[3] );
1098       gp_XYZ  v0 = p0t - p0s;
1099       gp_XYZ  v1 = p1t - p1s;
1100       gp_XYZ v01 = p1s - p0s;
1101       gp_XYZ   n = ( v0 ^ v01 ) + ( v1 ^ v01 );
1102       gp_XYZ   d = v01 ^ n;
1103       d.Normalize();
1104       return d;
1105     }
1106     gp_XYZ GetDir(_LayerEdge* le1, _LayerEdge* le2) // return average direction of _LayerEdge's
1107     {
1108       _nn[0]=le1->_nodes[0];
1109       _nn[1]=le1->_nodes.back();
1110       _nn[2]=le2->_nodes.back();
1111       _nn[3]=le2->_nodes[0];
1112       return GetDir();
1113     }
1114   };
1115   //--------------------------------------------------------------------------------
1116   /*!
1117    * \brief Retriever of node coordinates either directly or from a surface by node UV.
1118    * \warning Location of a surface is ignored
1119    */
1120   struct _NodeCoordHelper
1121   {
1122     SMESH_MesherHelper&        _helper;
1123     const TopoDS_Face&         _face;
1124     Handle(Geom_Surface)       _surface;
1125     gp_XYZ (_NodeCoordHelper::* _fun)(const SMDS_MeshNode* n) const;
1126
1127     _NodeCoordHelper(const TopoDS_Face& F, SMESH_MesherHelper& helper, bool is2D)
1128       : _helper( helper ), _face( F )
1129     {
1130       if ( is2D )
1131       {
1132         TopLoc_Location loc;
1133         _surface = BRep_Tool::Surface( _face, loc );
1134       }
1135       if ( _surface.IsNull() )
1136         _fun = & _NodeCoordHelper::direct;
1137       else
1138         _fun = & _NodeCoordHelper::byUV;
1139     }
1140     gp_XYZ operator()(const SMDS_MeshNode* n) const { return (this->*_fun)( n ); }
1141
1142   private:
1143     gp_XYZ direct(const SMDS_MeshNode* n) const
1144     {
1145       return SMESH_TNodeXYZ( n );
1146     }
1147     gp_XYZ byUV  (const SMDS_MeshNode* n) const
1148     {
1149       gp_XY uv = _helper.GetNodeUV( _face, n );
1150       return _surface->Value( uv.X(), uv.Y() ).XYZ();
1151     }
1152   };
1153
1154   //================================================================================
1155   /*!
1156    * \brief Check angle between vectors 
1157    */
1158   //================================================================================
1159
1160   inline bool isLessAngle( const gp_Vec& v1, const gp_Vec& v2, const double cos )
1161   {
1162     double dot = v1 * v2; // cos * |v1| * |v2|
1163     double l1  = v1.SquareMagnitude();
1164     double l2  = v2.SquareMagnitude();
1165     return (( dot * cos >= 0 ) && 
1166             ( dot * dot ) / l1 / l2 >= ( cos * cos ));
1167   }
1168
1169 } // namespace VISCOUS_3D
1170
1171
1172
1173 //================================================================================
1174 // StdMeshers_ViscousLayers hypothesis
1175 //
1176 StdMeshers_ViscousLayers::StdMeshers_ViscousLayers(int hypId, int studyId, SMESH_Gen* gen)
1177   :SMESH_Hypothesis(hypId, studyId, gen),
1178    _isToIgnoreShapes(1), _nbLayers(1), _thickness(1), _stretchFactor(1),
1179    _method( SURF_OFFSET_SMOOTH )
1180 {
1181   _name = StdMeshers_ViscousLayers::GetHypType();
1182   _param_algo_dim = -3; // auxiliary hyp used by 3D algos
1183 } // --------------------------------------------------------------------------------
1184 void StdMeshers_ViscousLayers::SetBndShapes(const std::vector<int>& faceIds, bool toIgnore)
1185 {
1186   if ( faceIds != _shapeIds )
1187     _shapeIds = faceIds, NotifySubMeshesHypothesisModification();
1188   if ( _isToIgnoreShapes != toIgnore )
1189     _isToIgnoreShapes = toIgnore, NotifySubMeshesHypothesisModification();
1190 } // --------------------------------------------------------------------------------
1191 void StdMeshers_ViscousLayers::SetTotalThickness(double thickness)
1192 {
1193   if ( thickness != _thickness )
1194     _thickness = thickness, NotifySubMeshesHypothesisModification();
1195 } // --------------------------------------------------------------------------------
1196 void StdMeshers_ViscousLayers::SetNumberLayers(int nb)
1197 {
1198   if ( _nbLayers != nb )
1199     _nbLayers = nb, NotifySubMeshesHypothesisModification();
1200 } // --------------------------------------------------------------------------------
1201 void StdMeshers_ViscousLayers::SetStretchFactor(double factor)
1202 {
1203   if ( _stretchFactor != factor )
1204     _stretchFactor = factor, NotifySubMeshesHypothesisModification();
1205 } // --------------------------------------------------------------------------------
1206 void StdMeshers_ViscousLayers::SetMethod( ExtrusionMethod method )
1207 {
1208   if ( _method != method )
1209     _method = method, NotifySubMeshesHypothesisModification();
1210 } // --------------------------------------------------------------------------------
1211 SMESH_ProxyMesh::Ptr
1212 StdMeshers_ViscousLayers::Compute(SMESH_Mesh&         theMesh,
1213                                   const TopoDS_Shape& theShape,
1214                                   const bool          toMakeN2NMap) const
1215 {
1216   using namespace VISCOUS_3D;
1217   _ViscousBuilder bulder;
1218   SMESH_ComputeErrorPtr err = bulder.Compute( theMesh, theShape );
1219   if ( err && !err->IsOK() )
1220     return SMESH_ProxyMesh::Ptr();
1221
1222   vector<SMESH_ProxyMesh::Ptr> components;
1223   TopExp_Explorer exp( theShape, TopAbs_SOLID );
1224   for ( ; exp.More(); exp.Next() )
1225   {
1226     if ( _MeshOfSolid* pm =
1227          _ViscousListener::GetSolidMesh( &theMesh, exp.Current(), /*toCreate=*/false))
1228     {
1229       if ( toMakeN2NMap && !pm->_n2nMapComputed )
1230         if ( !bulder.MakeN2NMap( pm ))
1231           return SMESH_ProxyMesh::Ptr();
1232       components.push_back( SMESH_ProxyMesh::Ptr( pm ));
1233       pm->myIsDeletable = false; // it will de deleted by boost::shared_ptr
1234
1235       if ( pm->_warning && !pm->_warning->IsOK() )
1236       {
1237         SMESH_subMesh* sm = theMesh.GetSubMesh( exp.Current() );
1238         SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1239         if ( !smError || smError->IsOK() )
1240           smError = pm->_warning;
1241       }
1242     }
1243     _ViscousListener::RemoveSolidMesh ( &theMesh, exp.Current() );
1244   }
1245   switch ( components.size() )
1246   {
1247   case 0: break;
1248
1249   case 1: return components[0];
1250
1251   default: return SMESH_ProxyMesh::Ptr( new SMESH_ProxyMesh( components ));
1252   }
1253   return SMESH_ProxyMesh::Ptr();
1254 } // --------------------------------------------------------------------------------
1255 std::ostream & StdMeshers_ViscousLayers::SaveTo(std::ostream & save)
1256 {
1257   save << " " << _nbLayers
1258        << " " << _thickness
1259        << " " << _stretchFactor
1260        << " " << _shapeIds.size();
1261   for ( size_t i = 0; i < _shapeIds.size(); ++i )
1262     save << " " << _shapeIds[i];
1263   save << " " << !_isToIgnoreShapes; // negate to keep the behavior in old studies.
1264   save << " " << _method;
1265   return save;
1266 } // --------------------------------------------------------------------------------
1267 std::istream & StdMeshers_ViscousLayers::LoadFrom(std::istream & load)
1268 {
1269   int nbFaces, faceID, shapeToTreat, method;
1270   load >> _nbLayers >> _thickness >> _stretchFactor >> nbFaces;
1271   while ( (int) _shapeIds.size() < nbFaces && load >> faceID )
1272     _shapeIds.push_back( faceID );
1273   if ( load >> shapeToTreat ) {
1274     _isToIgnoreShapes = !shapeToTreat;
1275     if ( load >> method )
1276       _method = (ExtrusionMethod) method;
1277   }
1278   else {
1279     _isToIgnoreShapes = true; // old behavior
1280   }
1281   return load;
1282 } // --------------------------------------------------------------------------------
1283 bool StdMeshers_ViscousLayers::SetParametersByMesh(const SMESH_Mesh*   theMesh,
1284                                                    const TopoDS_Shape& theShape)
1285 {
1286   // TODO
1287   return false;
1288 } // --------------------------------------------------------------------------------
1289 SMESH_ComputeErrorPtr
1290 StdMeshers_ViscousLayers::CheckHypothesis(SMESH_Mesh&                          theMesh,
1291                                           const TopoDS_Shape&                  theShape,
1292                                           SMESH_Hypothesis::Hypothesis_Status& theStatus)
1293 {
1294   VISCOUS_3D::_ViscousBuilder bulder;
1295   SMESH_ComputeErrorPtr err = bulder.CheckHypotheses( theMesh, theShape );
1296   if ( err && !err->IsOK() )
1297     theStatus = SMESH_Hypothesis::HYP_INCOMPAT_HYPS;
1298   else
1299     theStatus = SMESH_Hypothesis::HYP_OK;
1300
1301   return err;
1302 }
1303 // --------------------------------------------------------------------------------
1304 bool StdMeshers_ViscousLayers::IsShapeWithLayers(int shapeIndex) const
1305 {
1306   bool isIn =
1307     ( std::find( _shapeIds.begin(), _shapeIds.end(), shapeIndex ) != _shapeIds.end() );
1308   return IsToIgnoreShapes() ? !isIn : isIn;
1309 }
1310 // END StdMeshers_ViscousLayers hypothesis
1311 //================================================================================
1312
1313 namespace VISCOUS_3D
1314 {
1315   gp_XYZ getEdgeDir( const TopoDS_Edge& E, const TopoDS_Vertex& fromV )
1316   {
1317     gp_Vec dir;
1318     double f,l;
1319     Handle(Geom_Curve) c = BRep_Tool::Curve( E, f, l );
1320     if ( c.IsNull() ) return gp_XYZ( Precision::Infinite(), 1e100, 1e100 );
1321     gp_Pnt p = BRep_Tool::Pnt( fromV );
1322     double distF = p.SquareDistance( c->Value( f ));
1323     double distL = p.SquareDistance( c->Value( l ));
1324     c->D1(( distF < distL ? f : l), p, dir );
1325     if ( distL < distF ) dir.Reverse();
1326     return dir.XYZ();
1327   }
1328   //--------------------------------------------------------------------------------
1329   gp_XYZ getEdgeDir( const TopoDS_Edge& E, const SMDS_MeshNode* atNode,
1330                      SMESH_MesherHelper& helper)
1331   {
1332     gp_Vec dir;
1333     double f,l; gp_Pnt p;
1334     Handle(Geom_Curve) c = BRep_Tool::Curve( E, f, l );
1335     if ( c.IsNull() ) return gp_XYZ( Precision::Infinite(), 1e100, 1e100 );
1336     double u = helper.GetNodeU( E, atNode );
1337     c->D1( u, p, dir );
1338     return dir.XYZ();
1339   }
1340   //--------------------------------------------------------------------------------
1341   gp_XYZ getFaceDir( const TopoDS_Face& F, const TopoDS_Vertex& fromV,
1342                      const SMDS_MeshNode* node, SMESH_MesherHelper& helper, bool& ok,
1343                      double* cosin=0);
1344   //--------------------------------------------------------------------------------
1345   gp_XYZ getFaceDir( const TopoDS_Face& F, const TopoDS_Edge& fromE,
1346                      const SMDS_MeshNode* node, SMESH_MesherHelper& helper, bool& ok)
1347   {
1348     double f,l;
1349     Handle(Geom_Curve) c = BRep_Tool::Curve( fromE, f, l );
1350     if ( c.IsNull() )
1351     {
1352       TopoDS_Vertex v = helper.IthVertex( 0, fromE );
1353       return getFaceDir( F, v, node, helper, ok );
1354     }
1355     gp_XY uv = helper.GetNodeUV( F, node, 0, &ok );
1356     Handle(Geom_Surface) surface = BRep_Tool::Surface( F );
1357     gp_Pnt p; gp_Vec du, dv, norm;
1358     surface->D1( uv.X(),uv.Y(), p, du,dv );
1359     norm = du ^ dv;
1360
1361     double u = helper.GetNodeU( fromE, node, 0, &ok );
1362     c->D1( u, p, du );
1363     TopAbs_Orientation o = helper.GetSubShapeOri( F.Oriented(TopAbs_FORWARD), fromE);
1364     if ( o == TopAbs_REVERSED )
1365       du.Reverse();
1366
1367     gp_Vec dir = norm ^ du;
1368
1369     if ( node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX &&
1370          helper.IsClosedEdge( fromE ))
1371     {
1372       if ( fabs(u-f) < fabs(u-l)) c->D1( l, p, dv );
1373       else                        c->D1( f, p, dv );
1374       if ( o == TopAbs_REVERSED )
1375         dv.Reverse();
1376       gp_Vec dir2 = norm ^ dv;
1377       dir = dir.Normalized() + dir2.Normalized();
1378     }
1379     return dir.XYZ();
1380   }
1381   //--------------------------------------------------------------------------------
1382   gp_XYZ getFaceDir( const TopoDS_Face& F, const TopoDS_Vertex& fromV,
1383                      const SMDS_MeshNode* node, SMESH_MesherHelper& helper,
1384                      bool& ok, double* cosin)
1385   {
1386     TopoDS_Face faceFrw = F;
1387     faceFrw.Orientation( TopAbs_FORWARD );
1388     //double f,l; TopLoc_Location loc;
1389     TopoDS_Edge edges[2]; // sharing a vertex
1390     size_t nbEdges = 0;
1391     {
1392       TopoDS_Vertex VV[2];
1393       TopExp_Explorer exp( faceFrw, TopAbs_EDGE );
1394       for ( ; exp.More() && nbEdges < 2; exp.Next() )
1395       {
1396         const TopoDS_Edge& e = TopoDS::Edge( exp.Current() );
1397         if ( SMESH_Algo::isDegenerated( e )) continue;
1398         TopExp::Vertices( e, VV[0], VV[1], /*CumOri=*/true );
1399         if ( VV[1].IsSame( fromV )) {
1400           nbEdges += edges[ 0 ].IsNull();
1401           edges[ 0 ] = e;
1402         }
1403         else if ( VV[0].IsSame( fromV )) {
1404           nbEdges += edges[ 1 ].IsNull();
1405           edges[ 1 ] = e;
1406         }
1407       }
1408     }
1409     gp_XYZ dir(0,0,0), edgeDir[2];
1410     if ( nbEdges == 2 )
1411     {
1412       // get dirs of edges going fromV
1413       ok = true;
1414       for ( size_t i = 0; i < nbEdges && ok; ++i )
1415       {
1416         edgeDir[i] = getEdgeDir( edges[i], fromV );
1417         double size2 = edgeDir[i].SquareModulus();
1418         if (( ok = size2 > numeric_limits<double>::min() ))
1419           edgeDir[i] /= sqrt( size2 );
1420       }
1421       if ( !ok ) return dir;
1422
1423       // get angle between the 2 edges
1424       gp_Vec faceNormal;
1425       double angle = helper.GetAngle( edges[0], edges[1], faceFrw, fromV, &faceNormal );
1426       if ( Abs( angle ) < 5 * M_PI/180 )
1427       {
1428         dir = ( faceNormal.XYZ() ^ edgeDir[0].Reversed()) + ( faceNormal.XYZ() ^ edgeDir[1] );
1429       }
1430       else
1431       {
1432         dir = edgeDir[0] + edgeDir[1];
1433         if ( angle < 0 )
1434           dir.Reverse();
1435       }
1436       if ( cosin ) {
1437         double angle = gp_Vec( edgeDir[0] ).Angle( dir );
1438         *cosin = Cos( angle );
1439       }
1440     }
1441     else if ( nbEdges == 1 )
1442     {
1443       dir = getFaceDir( faceFrw, edges[ edges[0].IsNull() ], node, helper, ok );
1444       if ( cosin ) *cosin = 1.;
1445     }
1446     else
1447     {
1448       ok = false;
1449     }
1450
1451     return dir;
1452   }
1453
1454   //================================================================================
1455   /*!
1456    * \brief Finds concave VERTEXes of a FACE
1457    */
1458   //================================================================================
1459
1460   bool getConcaveVertices( const TopoDS_Face&  F,
1461                            SMESH_MesherHelper& helper,
1462                            set< TGeomID >*     vertices = 0)
1463   {
1464     // check angles at VERTEXes
1465     TError error;
1466     TSideVector wires = StdMeshers_FaceSide::GetFaceWires( F, *helper.GetMesh(), 0, error );
1467     for ( size_t iW = 0; iW < wires.size(); ++iW )
1468     {
1469       const int nbEdges = wires[iW]->NbEdges();
1470       if ( nbEdges < 2 && SMESH_Algo::isDegenerated( wires[iW]->Edge(0)))
1471         continue;
1472       for ( int iE1 = 0; iE1 < nbEdges; ++iE1 )
1473       {
1474         if ( SMESH_Algo::isDegenerated( wires[iW]->Edge( iE1 ))) continue;
1475         int iE2 = ( iE1 + 1 ) % nbEdges;
1476         while ( SMESH_Algo::isDegenerated( wires[iW]->Edge( iE2 )))
1477           iE2 = ( iE2 + 1 ) % nbEdges;
1478         TopoDS_Vertex V = wires[iW]->FirstVertex( iE2 );
1479         double angle = helper.GetAngle( wires[iW]->Edge( iE1 ),
1480                                         wires[iW]->Edge( iE2 ), F, V );
1481         if ( angle < -5. * M_PI / 180. )
1482         {
1483           if ( !vertices )
1484             return true;
1485           vertices->insert( helper.GetMeshDS()->ShapeToIndex( V ));
1486         }
1487       }
1488     }
1489     return vertices ? !vertices->empty() : false;
1490   }
1491
1492   //================================================================================
1493   /*!
1494    * \brief Returns true if a FACE is bound by a concave EDGE
1495    */
1496   //================================================================================
1497
1498   bool isConcave( const TopoDS_Face&  F,
1499                   SMESH_MesherHelper& helper,
1500                   set< TGeomID >*     vertices = 0 )
1501   {
1502     bool isConcv = false;
1503     // if ( helper.Count( F, TopAbs_WIRE, /*useMap=*/false) > 1 )
1504     //   return true;
1505     gp_Vec2d drv1, drv2;
1506     gp_Pnt2d p;
1507     TopExp_Explorer eExp( F.Oriented( TopAbs_FORWARD ), TopAbs_EDGE );
1508     for ( ; eExp.More(); eExp.Next() )
1509     {
1510       const TopoDS_Edge& E = TopoDS::Edge( eExp.Current() );
1511       if ( SMESH_Algo::isDegenerated( E )) continue;
1512       // check if 2D curve is concave
1513       BRepAdaptor_Curve2d curve( E, F );
1514       const int nbIntervals = curve.NbIntervals( GeomAbs_C2 );
1515       TColStd_Array1OfReal intervals(1, nbIntervals + 1 );
1516       curve.Intervals( intervals, GeomAbs_C2 );
1517       bool isConvex = true;
1518       for ( int i = 1; i <= nbIntervals && isConvex; ++i )
1519       {
1520         double u1 = intervals( i );
1521         double u2 = intervals( i+1 );
1522         curve.D2( 0.5*( u1+u2 ), p, drv1, drv2 );
1523         double cross = drv1 ^ drv2;
1524         if ( E.Orientation() == TopAbs_REVERSED )
1525           cross = -cross;
1526         isConvex = ( cross > -1e-9 ); // 0.1 );
1527       }
1528       if ( !isConvex )
1529       {
1530         //cout << "Concave FACE " << helper.GetMeshDS()->ShapeToIndex( F ) << endl;
1531         isConcv = true;
1532         if ( vertices )
1533           break;
1534         else
1535           return true;
1536       }
1537     }
1538
1539     // check angles at VERTEXes
1540     if ( getConcaveVertices( F, helper, vertices ))
1541       isConcv = true;
1542
1543     return isConcv;
1544   }
1545
1546   //================================================================================
1547   /*!
1548    * \brief Computes mimimal distance of face in-FACE nodes from an EDGE
1549    *  \param [in] face - the mesh face to treat
1550    *  \param [in] nodeOnEdge - a node on the EDGE
1551    *  \param [out] faceSize - the computed distance
1552    *  \return bool - true if faceSize computed
1553    */
1554   //================================================================================
1555
1556   bool getDistFromEdge( const SMDS_MeshElement* face,
1557                         const SMDS_MeshNode*    nodeOnEdge,
1558                         double &                faceSize )
1559   {
1560     faceSize = Precision::Infinite();
1561     bool done = false;
1562
1563     int nbN  = face->NbCornerNodes();
1564     int iOnE = face->GetNodeIndex( nodeOnEdge );
1565     int iNext[2] = { SMESH_MesherHelper::WrapIndex( iOnE+1, nbN ),
1566                      SMESH_MesherHelper::WrapIndex( iOnE-1, nbN ) };
1567     const SMDS_MeshNode* nNext[2] = { face->GetNode( iNext[0] ),
1568                                       face->GetNode( iNext[1] ) };
1569     gp_XYZ segVec, segEnd = SMESH_TNodeXYZ( nodeOnEdge ); // segment on EDGE
1570     double segLen = -1.;
1571     // look for two neighbor not in-FACE nodes of face
1572     for ( int i = 0; i < 2; ++i )
1573     {
1574       if ( nNext[i]->GetPosition()->GetDim() != 2 &&
1575            nNext[i]->GetID() < nodeOnEdge->GetID() )
1576       {
1577         // look for an in-FACE node
1578         for ( int iN = 0; iN < nbN; ++iN )
1579         {
1580           if ( iN == iOnE || iN == iNext[i] )
1581             continue;
1582           SMESH_TNodeXYZ pInFace = face->GetNode( iN );
1583           gp_XYZ v = pInFace - segEnd;
1584           if ( segLen < 0 )
1585           {
1586             segVec = SMESH_TNodeXYZ( nNext[i] ) - segEnd;
1587             segLen = segVec.Modulus();
1588           }
1589           double distToSeg = v.Crossed( segVec ).Modulus() / segLen;
1590           faceSize = Min( faceSize, distToSeg );
1591           done = true;
1592         }
1593         segLen = -1;
1594       }
1595     }
1596     return done;
1597   }
1598   //================================================================================
1599   /*!
1600    * \brief Return direction of axis or revolution of a surface
1601    */
1602   //================================================================================
1603
1604   bool getRovolutionAxis( const Adaptor3d_Surface& surface,
1605                           gp_Dir &                 axis )
1606   {
1607     switch ( surface.GetType() ) {
1608     case GeomAbs_Cone:
1609     {
1610       gp_Cone cone = surface.Cone();
1611       axis = cone.Axis().Direction();
1612       break;
1613     }
1614     case GeomAbs_Sphere:
1615     {
1616       gp_Sphere sphere = surface.Sphere();
1617       axis = sphere.Position().Direction();
1618       break;
1619     }
1620     case GeomAbs_SurfaceOfRevolution:
1621     {
1622       axis = surface.AxeOfRevolution().Direction();
1623       break;
1624     }
1625     //case GeomAbs_SurfaceOfExtrusion:
1626     case GeomAbs_OffsetSurface:
1627     {
1628       Handle(Adaptor3d_HSurface) base = surface.BasisSurface();
1629       return getRovolutionAxis( base->Surface(), axis );
1630     }
1631     default: return false;
1632     }
1633     return true;
1634   }
1635
1636   //--------------------------------------------------------------------------------
1637   // DEBUG. Dump intermediate node positions into a python script
1638   // HOWTO use: run python commands written in a console to see
1639   //  construction steps of viscous layers
1640 #ifdef __myDEBUG
1641   ofstream* py;
1642   int       theNbPyFunc;
1643   struct PyDump {
1644     PyDump(SMESH_Mesh& m) {
1645       int tag = 3 + m.GetId();
1646       const char* fname = "/tmp/viscous.py";
1647       cout << "execfile('"<<fname<<"')"<<endl;
1648       py = new ofstream(fname);
1649       *py << "import SMESH" << endl
1650           << "from salome.smesh import smeshBuilder" << endl
1651           << "smesh  = smeshBuilder.New(salome.myStudy)" << endl
1652           << "meshSO = smesh.GetCurrentStudy().FindObjectID('0:1:2:" << tag <<"')" << endl
1653           << "mesh   = smesh.Mesh( meshSO.GetObject() )"<<endl;
1654       theNbPyFunc = 0;
1655     }
1656     void Finish() {
1657       if (py) {
1658         *py << "mesh.GroupOnFilter(SMESH.VOLUME,'Viscous Prisms',"
1659           "smesh.GetFilter(SMESH.VOLUME,SMESH.FT_ElemGeomType,'=',SMESH.Geom_PENTA))"<<endl;
1660         *py << "mesh.GroupOnFilter(SMESH.VOLUME,'Neg Volumes',"
1661           "smesh.GetFilter(SMESH.VOLUME,SMESH.FT_Volume3D,'<',0))"<<endl;
1662       }
1663       delete py; py=0;
1664     }
1665     ~PyDump() { Finish(); cout << "NB FUNCTIONS: " << theNbPyFunc << endl; }
1666   };
1667 #define dumpFunction(f) { _dumpFunction(f, __LINE__);}
1668 #define dumpMove(n)     { _dumpMove(n, __LINE__);}
1669 #define dumpMoveComm(n,txt) { _dumpMove(n, __LINE__, txt);}
1670 #define dumpCmd(txt)    { _dumpCmd(txt, __LINE__);}
1671   void _dumpFunction(const string& fun, int ln)
1672   { if (py) *py<< "def "<<fun<<"(): # "<< ln <<endl; cout<<fun<<"()"<<endl; ++theNbPyFunc; }
1673   void _dumpMove(const SMDS_MeshNode* n, int ln, const char* txt="")
1674   { if (py) *py<< "  mesh.MoveNode( "<<n->GetID()<< ", "<< n->X()
1675                << ", "<<n->Y()<<", "<< n->Z()<< ")\t\t # "<< ln <<" "<< txt << endl; }
1676   void _dumpCmd(const string& txt, int ln)
1677   { if (py) *py<< "  "<<txt<<" # "<< ln <<endl; }
1678   void dumpFunctionEnd()
1679   { if (py) *py<< "  return"<< endl; }
1680   void dumpChangeNodes( const SMDS_MeshElement* f )
1681   { if (py) { *py<< "  mesh.ChangeElemNodes( " << f->GetID()<<", [";
1682       for ( int i=1; i < f->NbNodes(); ++i ) *py << f->GetNode(i-1)->GetID()<<", ";
1683       *py << f->GetNode( f->NbNodes()-1 )->GetID() << " ])"<< endl; }}
1684 #define debugMsg( txt ) { cout << "# "<< txt << " (line: " << __LINE__ << ")" << endl; }
1685
1686 #else
1687
1688   struct PyDump { PyDump(SMESH_Mesh&) {} void Finish() {} };
1689 #define dumpFunction(f) f
1690 #define dumpMove(n)
1691 #define dumpMoveComm(n,txt)
1692 #define dumpCmd(txt)
1693 #define dumpFunctionEnd()
1694 #define dumpChangeNodes(f) { if(f) {} } // prevent "unused variable 'f'" warning
1695 #define debugMsg( txt ) {}
1696
1697 #endif
1698 }
1699
1700 using namespace VISCOUS_3D;
1701
1702 //================================================================================
1703 /*!
1704  * \brief Constructor of _ViscousBuilder
1705  */
1706 //================================================================================
1707
1708 _ViscousBuilder::_ViscousBuilder()
1709 {
1710   _error = SMESH_ComputeError::New(COMPERR_OK);
1711   _tmpFaceID = 0;
1712 }
1713
1714 //================================================================================
1715 /*!
1716  * \brief Stores error description and returns false
1717  */
1718 //================================================================================
1719
1720 bool _ViscousBuilder::error(const string& text, int solidId )
1721 {
1722   const string prefix = string("Viscous layers builder: ");
1723   _error->myName    = COMPERR_ALGO_FAILED;
1724   _error->myComment = prefix + text;
1725   if ( _mesh )
1726   {
1727     SMESH_subMesh* sm = _mesh->GetSubMeshContaining( solidId );
1728     if ( !sm && !_sdVec.empty() )
1729       sm = _mesh->GetSubMeshContaining( solidId = _sdVec[0]._index );
1730     if ( sm && sm->GetSubShape().ShapeType() == TopAbs_SOLID )
1731     {
1732       SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1733       if ( smError && smError->myAlgo )
1734         _error->myAlgo = smError->myAlgo;
1735       smError = _error;
1736       sm->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
1737     }
1738     // set KO to all solids
1739     for ( size_t i = 0; i < _sdVec.size(); ++i )
1740     {
1741       if ( _sdVec[i]._index == solidId )
1742         continue;
1743       sm = _mesh->GetSubMesh( _sdVec[i]._solid );
1744       if ( !sm->IsEmpty() )
1745         continue;
1746       SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1747       if ( !smError || smError->IsOK() )
1748       {
1749         smError = SMESH_ComputeError::New( COMPERR_ALGO_FAILED, prefix + "failed");
1750         sm->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
1751       }
1752     }
1753   }
1754   makeGroupOfLE(); // debug
1755
1756   return false;
1757 }
1758
1759 //================================================================================
1760 /*!
1761  * \brief At study restoration, restore event listeners used to clear an inferior
1762  *  dim sub-mesh modified by viscous layers
1763  */
1764 //================================================================================
1765
1766 void _ViscousBuilder::RestoreListeners()
1767 {
1768   // TODO
1769 }
1770
1771 //================================================================================
1772 /*!
1773  * \brief computes SMESH_ProxyMesh::SubMesh::_n2n
1774  */
1775 //================================================================================
1776
1777 bool _ViscousBuilder::MakeN2NMap( _MeshOfSolid* pm )
1778 {
1779   SMESH_subMesh* solidSM = pm->mySubMeshes.front();
1780   TopExp_Explorer fExp( solidSM->GetSubShape(), TopAbs_FACE );
1781   for ( ; fExp.More(); fExp.Next() )
1782   {
1783     SMESHDS_SubMesh* srcSmDS = pm->GetMeshDS()->MeshElements( fExp.Current() );
1784     const SMESH_ProxyMesh::SubMesh* prxSmDS = pm->GetProxySubMesh( fExp.Current() );
1785
1786     if ( !srcSmDS || !prxSmDS || !srcSmDS->NbElements() || !prxSmDS->NbElements() )
1787       continue;
1788     if ( srcSmDS->GetElements()->next() == prxSmDS->GetElements()->next())
1789       continue;
1790
1791     if ( srcSmDS->NbElements() != prxSmDS->NbElements() )
1792       return error( "Different nb elements in a source and a proxy sub-mesh", solidSM->GetId());
1793
1794     SMDS_ElemIteratorPtr srcIt = srcSmDS->GetElements();
1795     SMDS_ElemIteratorPtr prxIt = prxSmDS->GetElements();
1796     while( prxIt->more() )
1797     {
1798       const SMDS_MeshElement* fSrc = srcIt->next();
1799       const SMDS_MeshElement* fPrx = prxIt->next();
1800       if ( fSrc->NbNodes() != fPrx->NbNodes())
1801         return error( "Different elements in a source and a proxy sub-mesh", solidSM->GetId());
1802       for ( int i = 0 ; i < fPrx->NbNodes(); ++i )
1803         pm->setNode2Node( fSrc->GetNode(i), fPrx->GetNode(i), prxSmDS );
1804     }
1805   }
1806   pm->_n2nMapComputed = true;
1807   return true;
1808 }
1809
1810 //================================================================================
1811 /*!
1812  * \brief Does its job
1813  */
1814 //================================================================================
1815
1816 SMESH_ComputeErrorPtr _ViscousBuilder::Compute(SMESH_Mesh&         theMesh,
1817                                                const TopoDS_Shape& theShape)
1818 {
1819   // TODO: set priority of solids during Gen::Compute()
1820
1821   _mesh = & theMesh;
1822
1823   // check if proxy mesh already computed
1824   TopExp_Explorer exp( theShape, TopAbs_SOLID );
1825   if ( !exp.More() )
1826     return error("No SOLID's in theShape"), _error;
1827
1828   if ( _ViscousListener::GetSolidMesh( _mesh, exp.Current(), /*toCreate=*/false))
1829     return SMESH_ComputeErrorPtr(); // everything already computed
1830
1831   PyDump debugDump( theMesh );
1832
1833   // TODO: ignore already computed SOLIDs 
1834   if ( !findSolidsWithLayers())
1835     return _error;
1836
1837   if ( !findFacesWithLayers() )
1838     return _error;
1839
1840   for ( size_t i = 0; i < _sdVec.size(); ++i )
1841   {
1842     if ( ! makeLayer(_sdVec[i]) )
1843       return _error;
1844
1845     if ( _sdVec[i]._n2eMap.size() == 0 )
1846       continue;
1847     
1848     if ( ! inflate(_sdVec[i]) )
1849       return _error;
1850
1851     if ( ! refine(_sdVec[i]) )
1852       return _error;
1853   }
1854   if ( !shrink() )
1855     return _error;
1856
1857   addBoundaryElements();
1858
1859   makeGroupOfLE(); // debug
1860   debugDump.Finish();
1861
1862   return _error;
1863 }
1864
1865 //================================================================================
1866 /*!
1867  * \brief Check validity of hypotheses
1868  */
1869 //================================================================================
1870
1871 SMESH_ComputeErrorPtr _ViscousBuilder::CheckHypotheses( SMESH_Mesh&         mesh,
1872                                                         const TopoDS_Shape& shape )
1873 {
1874   _mesh = & mesh;
1875
1876   if ( _ViscousListener::GetSolidMesh( _mesh, shape, /*toCreate=*/false))
1877     return SMESH_ComputeErrorPtr(); // everything already computed
1878
1879
1880   findSolidsWithLayers();
1881   bool ok = findFacesWithLayers( true );
1882
1883   // remove _MeshOfSolid's of _SolidData's
1884   for ( size_t i = 0; i < _sdVec.size(); ++i )
1885     _ViscousListener::RemoveSolidMesh( _mesh, _sdVec[i]._solid );
1886
1887   if ( !ok )
1888     return _error;
1889
1890   return SMESH_ComputeErrorPtr();
1891 }
1892
1893 //================================================================================
1894 /*!
1895  * \brief Finds SOLIDs to compute using viscous layers. Fills _sdVec
1896  */
1897 //================================================================================
1898
1899 bool _ViscousBuilder::findSolidsWithLayers()
1900 {
1901   // get all solids
1902   TopTools_IndexedMapOfShape allSolids;
1903   TopExp::MapShapes( _mesh->GetShapeToMesh(), TopAbs_SOLID, allSolids );
1904   _sdVec.reserve( allSolids.Extent());
1905
1906   SMESH_Gen* gen = _mesh->GetGen();
1907   SMESH_HypoFilter filter;
1908   for ( int i = 1; i <= allSolids.Extent(); ++i )
1909   {
1910     // find StdMeshers_ViscousLayers hyp assigned to the i-th solid
1911     SMESH_Algo* algo = gen->GetAlgo( *_mesh, allSolids(i) );
1912     if ( !algo ) continue;
1913     // TODO: check if algo is hidden
1914     const list <const SMESHDS_Hypothesis *> & allHyps =
1915       algo->GetUsedHypothesis(*_mesh, allSolids(i), /*ignoreAuxiliary=*/false);
1916     _SolidData* soData = 0;
1917     list< const SMESHDS_Hypothesis *>::const_iterator hyp = allHyps.begin();
1918     const StdMeshers_ViscousLayers* viscHyp = 0;
1919     for ( ; hyp != allHyps.end(); ++hyp )
1920       if (( viscHyp = dynamic_cast<const StdMeshers_ViscousLayers*>( *hyp )))
1921       {
1922         TopoDS_Shape hypShape;
1923         filter.Init( filter.Is( viscHyp ));
1924         _mesh->GetHypothesis( allSolids(i), filter, true, &hypShape );
1925
1926         if ( !soData )
1927         {
1928           _MeshOfSolid* proxyMesh = _ViscousListener::GetSolidMesh( _mesh,
1929                                                                     allSolids(i),
1930                                                                     /*toCreate=*/true);
1931           _sdVec.push_back( _SolidData( allSolids(i), proxyMesh ));
1932           soData = & _sdVec.back();
1933           soData->_index = getMeshDS()->ShapeToIndex( allSolids(i));
1934           soData->_helper = new SMESH_MesherHelper( *_mesh );
1935           soData->_helper->SetSubShape( allSolids(i) );
1936         }
1937         soData->_hyps.push_back( viscHyp );
1938         soData->_hypShapes.push_back( hypShape );
1939       }
1940   }
1941   if ( _sdVec.empty() )
1942     return error
1943       ( SMESH_Comment(StdMeshers_ViscousLayers::GetHypType()) << " hypothesis not found",0);
1944
1945   return true;
1946 }
1947
1948 //================================================================================
1949 /*!
1950  * \brief 
1951  */
1952 //================================================================================
1953
1954 bool _ViscousBuilder::findFacesWithLayers(const bool onlyWith)
1955 {
1956   SMESH_MesherHelper helper( *_mesh );
1957   TopExp_Explorer exp;
1958   TopTools_IndexedMapOfShape solids;
1959
1960   // collect all faces-to-ignore defined by hyp
1961   for ( size_t i = 0; i < _sdVec.size(); ++i )
1962   {
1963     solids.Add( _sdVec[i]._solid );
1964
1965     // get faces-to-ignore defined by each hyp
1966     typedef const StdMeshers_ViscousLayers* THyp;
1967     typedef std::pair< set<TGeomID>, THyp > TFacesOfHyp;
1968     list< TFacesOfHyp > ignoreFacesOfHyps;
1969     list< THyp >::iterator              hyp = _sdVec[i]._hyps.begin();
1970     list< TopoDS_Shape >::iterator hypShape = _sdVec[i]._hypShapes.begin();
1971     for ( ; hyp != _sdVec[i]._hyps.end(); ++hyp, ++hypShape )
1972     {
1973       ignoreFacesOfHyps.push_back( TFacesOfHyp( set<TGeomID>(), *hyp ));
1974       getIgnoreFaces( _sdVec[i]._solid, *hyp, *hypShape, ignoreFacesOfHyps.back().first );
1975     }
1976
1977     // fill _SolidData::_face2hyp and check compatibility of hypotheses
1978     const int nbHyps = _sdVec[i]._hyps.size();
1979     if ( nbHyps > 1 )
1980     {
1981       // check if two hypotheses define different parameters for the same FACE
1982       list< TFacesOfHyp >::iterator igFacesOfHyp;
1983       for ( exp.Init( _sdVec[i]._solid, TopAbs_FACE ); exp.More(); exp.Next() )
1984       {
1985         const TGeomID faceID = getMeshDS()->ShapeToIndex( exp.Current() );
1986         THyp hyp = 0;
1987         igFacesOfHyp = ignoreFacesOfHyps.begin();
1988         for ( ; igFacesOfHyp != ignoreFacesOfHyps.end(); ++igFacesOfHyp )
1989           if ( ! igFacesOfHyp->first.count( faceID ))
1990           {
1991             if ( hyp )
1992               return error(SMESH_Comment("Several hypotheses define "
1993                                          "Viscous Layers on the face #") << faceID );
1994             hyp = igFacesOfHyp->second;
1995           }
1996         if ( hyp )
1997           _sdVec[i]._face2hyp.insert( make_pair( faceID, hyp ));
1998         else
1999           _sdVec[i]._ignoreFaceIds.insert( faceID );
2000       }
2001
2002       // check if two hypotheses define different number of viscous layers for
2003       // adjacent faces of a solid
2004       set< int > nbLayersSet;
2005       igFacesOfHyp = ignoreFacesOfHyps.begin();
2006       for ( ; igFacesOfHyp != ignoreFacesOfHyps.end(); ++igFacesOfHyp )
2007       {
2008         nbLayersSet.insert( igFacesOfHyp->second->GetNumberLayers() );
2009       }
2010       if ( nbLayersSet.size() > 1 )
2011       {
2012         for ( exp.Init( _sdVec[i]._solid, TopAbs_EDGE ); exp.More(); exp.Next() )
2013         {
2014           PShapeIteratorPtr fIt = helper.GetAncestors( exp.Current(), *_mesh, TopAbs_FACE );
2015           THyp hyp1 = 0, hyp2 = 0;
2016           while( const TopoDS_Shape* face = fIt->next() )
2017           {
2018             const TGeomID faceID = getMeshDS()->ShapeToIndex( *face );
2019             map< TGeomID, THyp >::iterator f2h = _sdVec[i]._face2hyp.find( faceID );
2020             if ( f2h != _sdVec[i]._face2hyp.end() )
2021             {
2022               ( hyp1 ? hyp2 : hyp1 ) = f2h->second;
2023             }
2024           }
2025           if ( hyp1 && hyp2 &&
2026                hyp1->GetNumberLayers() != hyp2->GetNumberLayers() )
2027           {
2028             return error("Two hypotheses define different number of "
2029                          "viscous layers on adjacent faces");
2030           }
2031         }
2032       }
2033     } // if ( nbHyps > 1 )
2034     else
2035     {
2036       _sdVec[i]._ignoreFaceIds.swap( ignoreFacesOfHyps.back().first );
2037     }
2038   } // loop on _sdVec
2039
2040   if ( onlyWith ) // is called to check hypotheses compatibility only
2041     return true;
2042
2043   // fill _SolidData::_reversedFaceIds
2044   for ( size_t i = 0; i < _sdVec.size(); ++i )
2045   {
2046     exp.Init( _sdVec[i]._solid.Oriented( TopAbs_FORWARD ), TopAbs_FACE );
2047     for ( ; exp.More(); exp.Next() )
2048     {
2049       const TopoDS_Face& face = TopoDS::Face( exp.Current() );
2050       const TGeomID faceID = getMeshDS()->ShapeToIndex( face );
2051       if ( //!sdVec[i]._ignoreFaceIds.count( faceID ) &&
2052           helper.NbAncestors( face, *_mesh, TopAbs_SOLID ) > 1 &&
2053           helper.IsReversedSubMesh( face ))
2054       {
2055         _sdVec[i]._reversedFaceIds.insert( faceID );
2056       }
2057     }
2058   }
2059
2060   // Find faces to shrink mesh on (solution 2 in issue 0020832);
2061   TopTools_IndexedMapOfShape shapes;
2062   for ( size_t i = 0; i < _sdVec.size(); ++i )
2063   {
2064     shapes.Clear();
2065     TopExp::MapShapes(_sdVec[i]._solid, TopAbs_EDGE, shapes);
2066     for ( int iE = 1; iE <= shapes.Extent(); ++iE )
2067     {
2068       const TopoDS_Shape& edge = shapes(iE);
2069       // find 2 faces sharing an edge
2070       TopoDS_Shape FF[2];
2071       PShapeIteratorPtr fIt = helper.GetAncestors(edge, *_mesh, TopAbs_FACE);
2072       while ( fIt->more())
2073       {
2074         const TopoDS_Shape* f = fIt->next();
2075         if ( helper.IsSubShape( *f, _sdVec[i]._solid))
2076           FF[ int( !FF[0].IsNull()) ] = *f;
2077       }
2078       if( FF[1].IsNull() ) continue; // seam edge can be shared by 1 FACE only
2079       // check presence of layers on them
2080       int ignore[2];
2081       for ( int j = 0; j < 2; ++j )
2082         ignore[j] = _sdVec[i]._ignoreFaceIds.count ( getMeshDS()->ShapeToIndex( FF[j] ));
2083       if ( ignore[0] == ignore[1] )
2084         continue; // nothing interesting
2085       TopoDS_Shape fWOL = FF[ ignore[0] ? 0 : 1 ];
2086       // check presence of layers on fWOL within an adjacent SOLID
2087       bool collision = false;
2088       PShapeIteratorPtr sIt = helper.GetAncestors( fWOL, *_mesh, TopAbs_SOLID );
2089       while ( const TopoDS_Shape* solid = sIt->next() )
2090         if ( !solid->IsSame( _sdVec[i]._solid ))
2091         {
2092           int iSolid = solids.FindIndex( *solid );
2093           int  iFace = getMeshDS()->ShapeToIndex( fWOL );
2094           if ( iSolid > 0 && !_sdVec[ iSolid-1 ]._ignoreFaceIds.count( iFace ))
2095           {
2096             //_sdVec[i]._noShrinkShapes.insert( iFace );
2097             //fWOL.Nullify();
2098             collision = true;
2099           }
2100         }
2101       // add edge to maps
2102       if ( !fWOL.IsNull())
2103       {
2104         TGeomID edgeInd = getMeshDS()->ShapeToIndex( edge );
2105         _sdVec[i]._shrinkShape2Shape.insert( make_pair( edgeInd, fWOL ));
2106         if ( collision )
2107         {
2108           // _shrinkShape2Shape will be used to temporary inflate _LayerEdge's based
2109           // on the edge but shrink won't be performed
2110           _sdVec[i]._noShrinkShapes.insert( edgeInd );
2111         }
2112       }
2113     }
2114   }
2115   // Exclude from _shrinkShape2Shape FACE's that can't be shrinked since
2116   // the algo of the SOLID sharing the FACE does not support it
2117   set< string > notSupportAlgos; notSupportAlgos.insert("Hexa_3D");
2118   for ( size_t i = 0; i < _sdVec.size(); ++i )
2119   {
2120     map< TGeomID, TopoDS_Shape >::iterator e2f = _sdVec[i]._shrinkShape2Shape.begin();
2121     for ( ; e2f != _sdVec[i]._shrinkShape2Shape.end(); ++e2f )
2122     {
2123       const TopoDS_Shape& fWOL = e2f->second;
2124       const TGeomID     edgeID = e2f->first;
2125       bool notShrinkFace = false;
2126       PShapeIteratorPtr soIt = helper.GetAncestors(fWOL, *_mesh, TopAbs_SOLID);
2127       while ( soIt->more() )
2128       {
2129         const TopoDS_Shape* solid = soIt->next();
2130         if ( _sdVec[i]._solid.IsSame( *solid )) continue;
2131         SMESH_Algo* algo = _mesh->GetGen()->GetAlgo( *_mesh, *solid );
2132         if ( !algo || !notSupportAlgos.count( algo->GetName() )) continue;
2133         notShrinkFace = true;
2134         size_t iSolid = 0;
2135         for ( ; iSolid < _sdVec.size(); ++iSolid )
2136         {
2137           if ( _sdVec[iSolid]._solid.IsSame( *solid ) ) {
2138             if ( _sdVec[iSolid]._shrinkShape2Shape.count( edgeID ))
2139               notShrinkFace = false;
2140             break;
2141           }
2142         }
2143         if ( notShrinkFace )
2144         {
2145           _sdVec[i]._noShrinkShapes.insert( edgeID );
2146
2147           // add VERTEXes of the edge in _noShrinkShapes
2148           TopoDS_Shape edge = getMeshDS()->IndexToShape( edgeID );
2149           for ( TopoDS_Iterator vIt( edge ); vIt.More(); vIt.Next() )
2150             _sdVec[i]._noShrinkShapes.insert( getMeshDS()->ShapeToIndex( vIt.Value() ));
2151
2152           // check if there is a collision with to-shrink-from EDGEs in iSolid
2153           if ( iSolid == _sdVec.size() )
2154             continue; // no VL in the solid
2155           shapes.Clear();
2156           TopExp::MapShapes( fWOL, TopAbs_EDGE, shapes);
2157           for ( int iE = 1; iE <= shapes.Extent(); ++iE )
2158           {
2159             const TopoDS_Edge& E = TopoDS::Edge( shapes( iE ));
2160             const TGeomID    eID = getMeshDS()->ShapeToIndex( E );
2161             if ( eID == edgeID ||
2162                  !_sdVec[iSolid]._shrinkShape2Shape.count( eID ) ||
2163                  _sdVec[i]._noShrinkShapes.count( eID ))
2164               continue;
2165             for ( int is1st = 0; is1st < 2; ++is1st )
2166             {
2167               TopoDS_Vertex V = helper.IthVertex( is1st, E );
2168               if ( _sdVec[i]._noShrinkShapes.count( getMeshDS()->ShapeToIndex( V ) ))
2169               {
2170                 // _sdVec[i]._noShrinkShapes.insert( eID );
2171                 // V = helper.IthVertex( !is1st, E );
2172                 // _sdVec[i]._noShrinkShapes.insert( getMeshDS()->ShapeToIndex( V ));
2173                 //iE = 0; // re-start the loop on EDGEs of fWOL
2174                 return error("No way to make a conformal mesh with "
2175                              "the given set of faces with layers", _sdVec[i]._index);
2176               }
2177             }
2178           }
2179         }
2180
2181       } // while ( soIt->more() )
2182     } // loop on _sdVec[i]._shrinkShape2Shape
2183   } // loop on _sdVec to fill in _SolidData::_noShrinkShapes
2184
2185   // Find the SHAPE along which to inflate _LayerEdge based on VERTEX
2186
2187   for ( size_t i = 0; i < _sdVec.size(); ++i )
2188   {
2189     shapes.Clear();
2190     TopExp::MapShapes(_sdVec[i]._solid, TopAbs_VERTEX, shapes);
2191     for ( int iV = 1; iV <= shapes.Extent(); ++iV )
2192     {
2193       const TopoDS_Shape& vertex = shapes(iV);
2194       // find faces WOL sharing the vertex
2195       vector< TopoDS_Shape > facesWOL;
2196       size_t totalNbFaces = 0;
2197       PShapeIteratorPtr fIt = helper.GetAncestors(vertex, *_mesh, TopAbs_FACE);
2198       while ( fIt->more())
2199       {
2200         const TopoDS_Shape* f = fIt->next();
2201         if ( helper.IsSubShape( *f, _sdVec[i]._solid ) )
2202         {
2203           totalNbFaces++;
2204           const int fID = getMeshDS()->ShapeToIndex( *f );
2205           if ( _sdVec[i]._ignoreFaceIds.count ( fID ) /*&&
2206                !_sdVec[i]._noShrinkShapes.count( fID )*/)
2207             facesWOL.push_back( *f );
2208         }
2209       }
2210       if ( facesWOL.size() == totalNbFaces || facesWOL.empty() )
2211         continue; // no layers at this vertex or no WOL
2212       TGeomID vInd = getMeshDS()->ShapeToIndex( vertex );
2213       switch ( facesWOL.size() )
2214       {
2215       case 1:
2216       {
2217         helper.SetSubShape( facesWOL[0] );
2218         if ( helper.IsRealSeam( vInd )) // inflate along a seam edge?
2219         {
2220           TopoDS_Shape seamEdge;
2221           PShapeIteratorPtr eIt = helper.GetAncestors(vertex, *_mesh, TopAbs_EDGE);
2222           while ( eIt->more() && seamEdge.IsNull() )
2223           {
2224             const TopoDS_Shape* e = eIt->next();
2225             if ( helper.IsRealSeam( *e ) )
2226               seamEdge = *e;
2227           }
2228           if ( !seamEdge.IsNull() )
2229           {
2230             _sdVec[i]._shrinkShape2Shape.insert( make_pair( vInd, seamEdge ));
2231             break;
2232           }
2233         }
2234         _sdVec[i]._shrinkShape2Shape.insert( make_pair( vInd, facesWOL[0] ));
2235         break;
2236       }
2237       case 2:
2238       {
2239         // find an edge shared by 2 faces
2240         PShapeIteratorPtr eIt = helper.GetAncestors(vertex, *_mesh, TopAbs_EDGE);
2241         while ( eIt->more())
2242         {
2243           const TopoDS_Shape* e = eIt->next();
2244           if ( helper.IsSubShape( *e, facesWOL[0]) &&
2245                helper.IsSubShape( *e, facesWOL[1]))
2246           {
2247             _sdVec[i]._shrinkShape2Shape.insert( make_pair( vInd, *e )); break;
2248           }
2249         }
2250         break;
2251       }
2252       default:
2253         return error("Not yet supported case", _sdVec[i]._index);
2254       }
2255     }
2256   }
2257
2258   // add FACEs of other SOLIDs to _ignoreFaceIds
2259   for ( size_t i = 0; i < _sdVec.size(); ++i )
2260   {
2261     shapes.Clear();
2262     TopExp::MapShapes(_sdVec[i]._solid, TopAbs_FACE, shapes);
2263
2264     for ( exp.Init( _mesh->GetShapeToMesh(), TopAbs_FACE ); exp.More(); exp.Next() )
2265     {
2266       if ( !shapes.Contains( exp.Current() ))
2267         _sdVec[i]._ignoreFaceIds.insert( getMeshDS()->ShapeToIndex( exp.Current() ));
2268     }
2269   }
2270
2271   return true;
2272 }
2273
2274 //================================================================================
2275 /*!
2276  * \brief Finds FACEs w/o layers for a given SOLID by an hypothesis
2277  */
2278 //================================================================================
2279
2280 void _ViscousBuilder::getIgnoreFaces(const TopoDS_Shape&             solid,
2281                                      const StdMeshers_ViscousLayers* hyp,
2282                                      const TopoDS_Shape&             hypShape,
2283                                      set<TGeomID>&                   ignoreFaceIds)
2284 {
2285   TopExp_Explorer exp;
2286
2287   vector<TGeomID> ids = hyp->GetBndShapes();
2288   if ( hyp->IsToIgnoreShapes() ) // FACEs to ignore are given
2289   {
2290     for ( size_t ii = 0; ii < ids.size(); ++ii )
2291     {
2292       const TopoDS_Shape& s = getMeshDS()->IndexToShape( ids[ii] );
2293       if ( !s.IsNull() && s.ShapeType() == TopAbs_FACE )
2294         ignoreFaceIds.insert( ids[ii] );
2295     }
2296   }
2297   else // FACEs with layers are given
2298   {
2299     exp.Init( solid, TopAbs_FACE );
2300     for ( ; exp.More(); exp.Next() )
2301     {
2302       TGeomID faceInd = getMeshDS()->ShapeToIndex( exp.Current() );
2303       if ( find( ids.begin(), ids.end(), faceInd ) == ids.end() )
2304         ignoreFaceIds.insert( faceInd );
2305     }
2306   }
2307
2308   // ignore internal FACEs if inlets and outlets are specified
2309   if ( hyp->IsToIgnoreShapes() )
2310   {
2311     TopTools_IndexedDataMapOfShapeListOfShape solidsOfFace;
2312     TopExp::MapShapesAndAncestors( hypShape,
2313                                    TopAbs_FACE, TopAbs_SOLID, solidsOfFace);
2314
2315     for ( exp.Init( solid, TopAbs_FACE ); exp.More(); exp.Next() )
2316     {
2317       const TopoDS_Face& face = TopoDS::Face( exp.Current() );
2318       if ( SMESH_MesherHelper::NbAncestors( face, *_mesh, TopAbs_SOLID ) < 2 )
2319         continue;
2320
2321       int nbSolids = solidsOfFace.FindFromKey( face ).Extent();
2322       if ( nbSolids > 1 )
2323         ignoreFaceIds.insert( getMeshDS()->ShapeToIndex( face ));
2324     }
2325   }
2326 }
2327
2328 //================================================================================
2329 /*!
2330  * \brief Create the inner surface of the viscous layer and prepare data for infation
2331  */
2332 //================================================================================
2333
2334 bool _ViscousBuilder::makeLayer(_SolidData& data)
2335 {
2336   // get all sub-shapes to make layers on
2337   set<TGeomID> subIds, faceIds;
2338   subIds = data._noShrinkShapes;
2339   TopExp_Explorer exp( data._solid, TopAbs_FACE );
2340   for ( ; exp.More(); exp.Next() )
2341   {
2342     SMESH_subMesh* fSubM = _mesh->GetSubMesh( exp.Current() );
2343     if ( ! data._ignoreFaceIds.count( fSubM->GetId() ))
2344       faceIds.insert( fSubM->GetId() );
2345   }
2346
2347   // make a map to find new nodes on sub-shapes shared with other SOLID
2348   map< TGeomID, TNode2Edge* >::iterator s2ne;
2349   map< TGeomID, TopoDS_Shape >::iterator s2s = data._shrinkShape2Shape.begin();
2350   for (; s2s != data._shrinkShape2Shape.end(); ++s2s )
2351   {
2352     TGeomID shapeInd = s2s->first;
2353     for ( size_t i = 0; i < _sdVec.size(); ++i )
2354     {
2355       if ( _sdVec[i]._index == data._index ) continue;
2356       map< TGeomID, TopoDS_Shape >::iterator s2s2 = _sdVec[i]._shrinkShape2Shape.find( shapeInd );
2357       if ( s2s2 != _sdVec[i]._shrinkShape2Shape.end() &&
2358            *s2s == *s2s2 && !_sdVec[i]._n2eMap.empty() )
2359       {
2360         data._s2neMap.insert( make_pair( shapeInd, &_sdVec[i]._n2eMap ));
2361         break;
2362       }
2363     }
2364   }
2365
2366   // Create temporary faces and _LayerEdge's
2367
2368   dumpFunction(SMESH_Comment("makeLayers_")<<data._index);
2369
2370   data._stepSize = Precision::Infinite();
2371   data._stepSizeNodes[0] = 0;
2372
2373   SMESH_MesherHelper helper( *_mesh );
2374   helper.SetSubShape( data._solid );
2375   helper.SetElementsOnShape( true );
2376
2377   vector< const SMDS_MeshNode*> newNodes; // of a mesh face
2378   TNode2Edge::iterator n2e2;
2379
2380   // collect _LayerEdge's of shapes they are based on
2381   vector< _EdgesOnShape >& edgesByGeom = data._edgesOnShape;
2382   const int nbShapes = getMeshDS()->MaxShapeIndex();
2383   edgesByGeom.resize( nbShapes+1 );
2384
2385   // set data of _EdgesOnShape's
2386   if ( SMESH_subMesh* sm = _mesh->GetSubMesh( data._solid ))
2387   {
2388     SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(/*includeSelf=*/false);
2389     while ( smIt->more() )
2390     {
2391       sm = smIt->next();
2392       if ( sm->GetSubShape().ShapeType() == TopAbs_FACE &&
2393            !faceIds.count( sm->GetId() ))
2394         continue;
2395       setShapeData( edgesByGeom[ sm->GetId() ], sm, data );
2396     }
2397   }
2398   // make _LayerEdge's
2399   for ( set<TGeomID>::iterator id = faceIds.begin(); id != faceIds.end(); ++id )
2400   {
2401     const TopoDS_Face& F = TopoDS::Face( getMeshDS()->IndexToShape( *id ));
2402     SMESH_subMesh* sm = _mesh->GetSubMesh( F );
2403     SMESH_ProxyMesh::SubMesh* proxySub =
2404       data._proxyMesh->getFaceSubM( F, /*create=*/true);
2405
2406     SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
2407     if ( !smDS ) return error(SMESH_Comment("Not meshed face ") << *id, data._index );
2408
2409     SMDS_ElemIteratorPtr eIt = smDS->GetElements();
2410     while ( eIt->more() )
2411     {
2412       const SMDS_MeshElement* face = eIt->next();
2413       double          faceMaxCosin = -1;
2414       _LayerEdge*     maxCosinEdge = 0;
2415       int             nbDegenNodes = 0;
2416
2417       newNodes.resize( face->NbCornerNodes() );
2418       for ( size_t i = 0 ; i < newNodes.size(); ++i )
2419       {
2420         const SMDS_MeshNode* n = face->GetNode( i );
2421         const int      shapeID = n->getshapeId();
2422         const bool onDegenShap = helper.IsDegenShape( shapeID );
2423         const bool onDegenEdge = ( onDegenShap && n->GetPosition()->GetDim() == 1 );
2424         if ( onDegenShap )
2425         {
2426           if ( onDegenEdge )
2427           {
2428             // substitute n on a degenerated EDGE with a node on a corresponding VERTEX
2429             const TopoDS_Shape& E = getMeshDS()->IndexToShape( shapeID );
2430             TopoDS_Vertex       V = helper.IthVertex( 0, TopoDS::Edge( E ));
2431             if ( const SMDS_MeshNode* vN = SMESH_Algo::VertexNode( V, getMeshDS() )) {
2432               n = vN;
2433               nbDegenNodes++;
2434             }
2435           }
2436           else
2437           {
2438             nbDegenNodes++;
2439           }
2440         }
2441         TNode2Edge::iterator n2e = data._n2eMap.insert( make_pair( n, (_LayerEdge*)0 )).first;
2442         if ( !(*n2e).second )
2443         {
2444           // add a _LayerEdge
2445           _LayerEdge* edge = new _LayerEdge();
2446           edge->_nodes.push_back( n );
2447           n2e->second = edge;
2448           edgesByGeom[ shapeID ]._edges.push_back( edge );
2449           const bool noShrink = data._noShrinkShapes.count( shapeID );
2450
2451           SMESH_TNodeXYZ xyz( n );
2452
2453           // set edge data or find already refined _LayerEdge and get data from it
2454           if (( !noShrink                                                     ) &&
2455               ( n->GetPosition()->GetTypeOfPosition() != SMDS_TOP_FACE        ) &&
2456               (( s2ne = data._s2neMap.find( shapeID )) != data._s2neMap.end() ) &&
2457               (( n2e2 = (*s2ne).second->find( n )) != s2ne->second->end()     ))
2458           {
2459             _LayerEdge* foundEdge = (*n2e2).second;
2460             gp_XYZ        lastPos = edge->Copy( *foundEdge, edgesByGeom[ shapeID ], helper );
2461             foundEdge->_pos.push_back( lastPos );
2462             // location of the last node is modified and we restore it by foundEdge->_pos.back()
2463             const_cast< SMDS_MeshNode* >
2464               ( edge->_nodes.back() )->setXYZ( xyz.X(), xyz.Y(), xyz.Z() );
2465           }
2466           else
2467           {
2468             if ( !noShrink )
2469             {
2470               edge->_nodes.push_back( helper.AddNode( xyz.X(), xyz.Y(), xyz.Z() ));
2471             }
2472             if ( !setEdgeData( *edge, edgesByGeom[ shapeID ], helper, data ))
2473               return false;
2474
2475             if ( edge->_nodes.size() < 2 )
2476               edge->Block( data );
2477               //data._noShrinkShapes.insert( shapeID );
2478           }
2479           dumpMove(edge->_nodes.back());
2480
2481           if ( edge->_cosin > faceMaxCosin )
2482           {
2483             faceMaxCosin = edge->_cosin;
2484             maxCosinEdge = edge;
2485           }
2486         }
2487         newNodes[ i ] = n2e->second->_nodes.back();
2488
2489         if ( onDegenEdge )
2490           data._n2eMap.insert( make_pair( face->GetNode( i ), n2e->second ));
2491       }
2492       if ( newNodes.size() - nbDegenNodes < 2 )
2493         continue;
2494
2495       // create a temporary face
2496       const SMDS_MeshElement* newFace =
2497         new _TmpMeshFace( newNodes, --_tmpFaceID, face->getshapeId(), face->getIdInShape() );
2498       proxySub->AddElement( newFace );
2499
2500       // compute inflation step size by min size of element on a convex surface
2501       if ( faceMaxCosin > theMinSmoothCosin )
2502         limitStepSize( data, face, maxCosinEdge );
2503
2504     } // loop on 2D elements on a FACE
2505   } // loop on FACEs of a SOLID to create _LayerEdge's
2506
2507
2508   // Set _LayerEdge::_neibors
2509   TNode2Edge::iterator n2e;
2510   for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
2511   {
2512     _EdgesOnShape& eos = data._edgesOnShape[iS];
2513     for ( size_t i = 0; i < eos._edges.size(); ++i )
2514     {
2515       _LayerEdge* edge = eos._edges[i];
2516       TIDSortedNodeSet nearNodes;
2517       SMDS_ElemIteratorPtr fIt = edge->_nodes[0]->GetInverseElementIterator(SMDSAbs_Face);
2518       while ( fIt->more() )
2519       {
2520         const SMDS_MeshElement* f = fIt->next();
2521         if ( !data._ignoreFaceIds.count( f->getshapeId() ))
2522           nearNodes.insert( f->begin_nodes(), f->end_nodes() );
2523       }
2524       nearNodes.erase( edge->_nodes[0] );
2525       edge->_neibors.reserve( nearNodes.size() );
2526       TIDSortedNodeSet::iterator node = nearNodes.begin();
2527       for ( ; node != nearNodes.end(); ++node )
2528         if (( n2e = data._n2eMap.find( *node )) != data._n2eMap.end() )
2529           edge->_neibors.push_back( n2e->second );
2530     }
2531   }
2532
2533   data._epsilon = 1e-7;
2534   if ( data._stepSize < 1. )
2535     data._epsilon *= data._stepSize;
2536
2537   if ( !findShapesToSmooth( data ))
2538     return false;
2539
2540   // limit data._stepSize depending on surface curvature and fill data._convexFaces
2541   limitStepSizeByCurvature( data ); // !!! it must be before node substitution in _Simplex
2542
2543   // Set target nodes into _Simplex and _LayerEdge's to _2NearEdges
2544   const SMDS_MeshNode* nn[2];
2545   for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
2546   {
2547     _EdgesOnShape& eos = data._edgesOnShape[iS];
2548     for ( size_t i = 0; i < eos._edges.size(); ++i )
2549     {
2550       _LayerEdge* edge = eos._edges[i];
2551       if ( edge->IsOnEdge() )
2552       {
2553         // get neighbor nodes
2554         bool hasData = ( edge->_2neibors->_edges[0] );
2555         if ( hasData ) // _LayerEdge is a copy of another one
2556         {
2557           nn[0] = edge->_2neibors->srcNode(0);
2558           nn[1] = edge->_2neibors->srcNode(1);
2559         }
2560         else if ( !findNeiborsOnEdge( edge, nn[0],nn[1], eos, data ))
2561         {
2562           return false;
2563         }
2564         // set neighbor _LayerEdge's
2565         for ( int j = 0; j < 2; ++j )
2566         {
2567           if (( n2e = data._n2eMap.find( nn[j] )) == data._n2eMap.end() )
2568             return error("_LayerEdge not found by src node", data._index);
2569           edge->_2neibors->_edges[j] = n2e->second;
2570         }
2571         if ( !hasData )
2572           edge->SetDataByNeighbors( nn[0], nn[1], eos, helper );
2573       }
2574
2575       for ( size_t j = 0; j < edge->_simplices.size(); ++j )
2576       {
2577         _Simplex& s = edge->_simplices[j];
2578         s._nNext = data._n2eMap[ s._nNext ]->_nodes.back();
2579         s._nPrev = data._n2eMap[ s._nPrev ]->_nodes.back();
2580       }
2581
2582       // For an _LayerEdge on a degenerated EDGE, copy some data from
2583       // a corresponding _LayerEdge on a VERTEX
2584       // (issue 52453, pb on a downloaded SampleCase2-Tet-netgen-mephisto.hdf)
2585       if ( helper.IsDegenShape( edge->_nodes[0]->getshapeId() ))
2586       {
2587         // Generally we should not get here
2588         if ( eos.ShapeType() != TopAbs_EDGE )
2589           continue;
2590         TopoDS_Vertex V = helper.IthVertex( 0, TopoDS::Edge( eos._shape ));
2591         const SMDS_MeshNode* vN = SMESH_Algo::VertexNode( V, getMeshDS() );
2592         if (( n2e = data._n2eMap.find( vN )) == data._n2eMap.end() )
2593           continue;
2594         const _LayerEdge* vEdge = n2e->second;
2595         edge->_normal    = vEdge->_normal;
2596         edge->_lenFactor = vEdge->_lenFactor;
2597         edge->_cosin     = vEdge->_cosin;
2598       }
2599
2600     } // loop on data._edgesOnShape._edges
2601   } // loop on data._edgesOnShape
2602
2603   // fix _LayerEdge::_2neibors on EDGEs to smooth
2604   // map< TGeomID,Handle(Geom_Curve)>::iterator e2c = data._edge2curve.begin();
2605   // for ( ; e2c != data._edge2curve.end(); ++e2c )
2606   //   if ( !e2c->second.IsNull() )
2607   //   {
2608   //     if ( _EdgesOnShape* eos = data.GetShapeEdges( e2c->first ))
2609   //       data.Sort2NeiborsOnEdge( eos->_edges );
2610   //   }
2611
2612   dumpFunctionEnd();
2613   return true;
2614 }
2615
2616 //================================================================================
2617 /*!
2618  * \brief Compute inflation step size by min size of element on a convex surface
2619  */
2620 //================================================================================
2621
2622 void _ViscousBuilder::limitStepSize( _SolidData&             data,
2623                                      const SMDS_MeshElement* face,
2624                                      const _LayerEdge*       maxCosinEdge )
2625 {
2626   int iN = 0;
2627   double minSize = 10 * data._stepSize;
2628   const int nbNodes = face->NbCornerNodes();
2629   for ( int i = 0; i < nbNodes; ++i )
2630   {
2631     const SMDS_MeshNode* nextN = face->GetNode( SMESH_MesherHelper::WrapIndex( i+1, nbNodes ));
2632     const SMDS_MeshNode*  curN = face->GetNode( i );
2633     if ( nextN->GetPosition()->GetTypeOfPosition() == SMDS_TOP_FACE ||
2634          curN-> GetPosition()->GetTypeOfPosition() == SMDS_TOP_FACE )
2635     {
2636       double dist = SMESH_TNodeXYZ( curN ).Distance( nextN );
2637       if ( dist < minSize )
2638         minSize = dist, iN = i;
2639     }
2640   }
2641   double newStep = 0.8 * minSize / maxCosinEdge->_lenFactor;
2642   if ( newStep < data._stepSize )
2643   {
2644     data._stepSize = newStep;
2645     data._stepSizeCoeff = 0.8 / maxCosinEdge->_lenFactor;
2646     data._stepSizeNodes[0] = face->GetNode( iN );
2647     data._stepSizeNodes[1] = face->GetNode( SMESH_MesherHelper::WrapIndex( iN+1, nbNodes ));
2648   }
2649 }
2650
2651 //================================================================================
2652 /*!
2653  * \brief Compute inflation step size by min size of element on a convex surface
2654  */
2655 //================================================================================
2656
2657 void _ViscousBuilder::limitStepSize( _SolidData& data, const double minSize )
2658 {
2659   if ( minSize < data._stepSize )
2660   {
2661     data._stepSize = minSize;
2662     if ( data._stepSizeNodes[0] )
2663     {
2664       double dist =
2665         SMESH_TNodeXYZ(data._stepSizeNodes[0]).Distance(data._stepSizeNodes[1]);
2666       data._stepSizeCoeff = data._stepSize / dist;
2667     }
2668   }
2669 }
2670
2671 //================================================================================
2672 /*!
2673  * \brief Limit data._stepSize by evaluating curvature of shapes and fill data._convexFaces
2674  */
2675 //================================================================================
2676
2677 void _ViscousBuilder::limitStepSizeByCurvature( _SolidData& data )
2678 {
2679   const int nbTestPnt = 5; // on a FACE sub-shape
2680
2681   BRepLProp_SLProps surfProp( 2, 1e-6 );
2682   SMESH_MesherHelper helper( *_mesh );
2683
2684   data._convexFaces.clear();
2685
2686   for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
2687   {
2688     _EdgesOnShape& eof = data._edgesOnShape[iS];
2689     if ( eof.ShapeType() != TopAbs_FACE ||
2690          data._ignoreFaceIds.count( eof._shapeID ))
2691       continue;
2692
2693     TopoDS_Face        F = TopoDS::Face( eof._shape );
2694     SMESH_subMesh *   sm = eof._subMesh;
2695     const TGeomID faceID = eof._shapeID;
2696
2697     BRepAdaptor_Surface surface( F, false );
2698     surfProp.SetSurface( surface );
2699
2700     bool isTooCurved = false;
2701
2702     _ConvexFace cnvFace;
2703     const double        oriFactor = ( F.Orientation() == TopAbs_REVERSED ? +1. : -1. );
2704     SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(/*includeSelf=*/true);
2705     while ( smIt->more() )
2706     {
2707       sm = smIt->next();
2708       const TGeomID subID = sm->GetId();
2709       // find _LayerEdge's of a sub-shape
2710       _EdgesOnShape* eos;
2711       if (( eos = data.GetShapeEdges( subID )))
2712         cnvFace._subIdToEOS.insert( make_pair( subID, eos ));
2713       else
2714         continue;
2715       // check concavity and curvature and limit data._stepSize
2716       const double minCurvature =
2717         1. / ( eos->_hyp.GetTotalThickness() * ( 1 + theThickToIntersection ));
2718       size_t iStep = Max( 1, eos->_edges.size() / nbTestPnt );
2719       for ( size_t i = 0; i < eos->_edges.size(); i += iStep )
2720       {
2721         gp_XY uv = helper.GetNodeUV( F, eos->_edges[ i ]->_nodes[0] );
2722         surfProp.SetParameters( uv.X(), uv.Y() );
2723         if ( !surfProp.IsCurvatureDefined() )
2724           continue;
2725         if ( surfProp.MaxCurvature() * oriFactor > minCurvature )
2726         {
2727           limitStepSize( data, 0.9 / surfProp.MaxCurvature() * oriFactor );
2728           isTooCurved = true;
2729         }
2730         if ( surfProp.MinCurvature() * oriFactor > minCurvature )
2731         {
2732           limitStepSize( data, 0.9 / surfProp.MinCurvature() * oriFactor );
2733           isTooCurved = true;
2734         }
2735       }
2736     } // loop on sub-shapes of the FACE
2737
2738     if ( !isTooCurved ) continue;
2739
2740     _ConvexFace & convFace =
2741       data._convexFaces.insert( make_pair( faceID, cnvFace )).first->second;
2742
2743     convFace._face = F;
2744     convFace._normalsFixed = false;
2745
2746     // Fill _ConvexFace::_simplexTestEdges. These _LayerEdge's are used to detect
2747     // prism distortion.
2748     map< TGeomID, _EdgesOnShape* >::iterator id2eos = convFace._subIdToEOS.find( faceID );
2749     if ( id2eos != convFace._subIdToEOS.end() && !id2eos->second->_edges.empty() )
2750     {
2751       // there are _LayerEdge's on the FACE it-self;
2752       // select _LayerEdge's near EDGEs
2753       _EdgesOnShape& eos = * id2eos->second;
2754       for ( size_t i = 0; i < eos._edges.size(); ++i )
2755       {
2756         _LayerEdge* ledge = eos._edges[ i ];
2757         for ( size_t j = 0; j < ledge->_simplices.size(); ++j )
2758           if ( ledge->_simplices[j]._nNext->GetPosition()->GetDim() < 2 )
2759           {
2760             convFace._simplexTestEdges.push_back( ledge );
2761             break;
2762           }
2763       }
2764     }
2765     else
2766     {
2767       // where there are no _LayerEdge's on a _ConvexFace,
2768       // as e.g. on a fillet surface with no internal nodes - issue 22580,
2769       // so that collision of viscous internal faces is not detected by check of
2770       // intersection of _LayerEdge's with the viscous internal faces.
2771
2772       set< const SMDS_MeshNode* > usedNodes;
2773
2774       // look for _LayerEdge's with null _sWOL
2775       id2eos = convFace._subIdToEOS.begin();
2776       for ( ; id2eos != convFace._subIdToEOS.end(); ++id2eos )
2777       {
2778         _EdgesOnShape& eos = * id2eos->second;
2779         if ( !eos._sWOL.IsNull() )
2780           continue;
2781         for ( size_t i = 0; i < eos._edges.size(); ++i )
2782         {
2783           _LayerEdge* ledge = eos._edges[ i ];
2784           const SMDS_MeshNode* srcNode = ledge->_nodes[0];
2785           if ( !usedNodes.insert( srcNode ).second ) continue;
2786
2787           for ( size_t i = 0; i < ledge->_simplices.size(); ++i )
2788           {
2789             usedNodes.insert( ledge->_simplices[i]._nPrev );
2790             usedNodes.insert( ledge->_simplices[i]._nNext );
2791           }
2792           convFace._simplexTestEdges.push_back( ledge );
2793         }
2794       }
2795     }
2796   } // loop on FACEs of data._solid
2797 }
2798
2799 //================================================================================
2800 /*!
2801  * \brief Detect shapes (and _LayerEdge's on them) to smooth
2802  */
2803 //================================================================================
2804
2805 bool _ViscousBuilder::findShapesToSmooth( _SolidData& data )
2806 {
2807   // define allowed thickness
2808   computeGeomSize( data ); // compute data._geomSize and _LayerEdge::_maxLen
2809
2810   data._maxThickness = 0;
2811   data._minThickness = 1e100;
2812   list< const StdMeshers_ViscousLayers* >::iterator hyp = data._hyps.begin();
2813   for ( ; hyp != data._hyps.end(); ++hyp )
2814   {
2815     data._maxThickness = Max( data._maxThickness, (*hyp)->GetTotalThickness() );
2816     data._minThickness = Min( data._minThickness, (*hyp)->GetTotalThickness() );
2817   }
2818   //const double tgtThick = /*Min( 0.5 * data._geomSize, */data._maxThickness;
2819
2820   // Find shapes needing smoothing; such a shape has _LayerEdge._normal on it's
2821   // boundry inclined to the shape at a sharp angle
2822
2823   //list< TGeomID > shapesToSmooth;
2824   TopTools_MapOfShape edgesOfSmooFaces;
2825
2826   SMESH_MesherHelper helper( *_mesh );
2827   bool ok = true;
2828
2829   vector< _EdgesOnShape >& edgesByGeom = data._edgesOnShape;
2830   data._nbShapesToSmooth = 0;
2831
2832   for ( size_t iS = 0; iS < edgesByGeom.size(); ++iS ) // check FACEs
2833   {
2834     _EdgesOnShape& eos = edgesByGeom[iS];
2835     eos._toSmooth = false;
2836     if ( eos._edges.empty() || eos.ShapeType() != TopAbs_FACE )
2837       continue;
2838
2839     double tgtThick = eos._hyp.GetTotalThickness();
2840     TopExp_Explorer eExp( edgesByGeom[iS]._shape, TopAbs_EDGE );
2841     for ( ; eExp.More() && !eos._toSmooth; eExp.Next() )
2842     {
2843       TGeomID iE = getMeshDS()->ShapeToIndex( eExp.Current() );
2844       vector<_LayerEdge*>& eE = edgesByGeom[ iE ]._edges;
2845       if ( eE.empty() ) continue;
2846
2847       double faceSize;
2848       for ( size_t i = 0; i < eE.size() && !eos._toSmooth; ++i )
2849         if ( eE[i]->_cosin > theMinSmoothCosin )
2850         {
2851           SMDS_ElemIteratorPtr fIt = eE[i]->_nodes[0]->GetInverseElementIterator(SMDSAbs_Face);
2852           while ( fIt->more() && !eos._toSmooth )
2853           {
2854             const SMDS_MeshElement* face = fIt->next();
2855             if ( face->getshapeId() == eos._shapeID &&
2856                  getDistFromEdge( face, eE[i]->_nodes[0], faceSize ))
2857             {
2858               eos._toSmooth = needSmoothing( eE[i]->_cosin, tgtThick, faceSize );
2859             }
2860           }
2861         }
2862     }
2863     if ( eos._toSmooth )
2864     {
2865       for ( eExp.ReInit(); eExp.More(); eExp.Next() )
2866         edgesOfSmooFaces.Add( eExp.Current() );
2867
2868       data.PrepareEdgesToSmoothOnFace( &edgesByGeom[iS], /*substituteSrcNodes=*/false );
2869     }
2870     data._nbShapesToSmooth += eos._toSmooth;
2871
2872   }  // check FACEs
2873
2874   for ( size_t iS = 0; iS < edgesByGeom.size(); ++iS ) // check EDGEs
2875   {
2876     _EdgesOnShape& eos = edgesByGeom[iS];
2877     eos._edgeSmoother = NULL;
2878     if ( eos._edges.empty() || eos.ShapeType() != TopAbs_EDGE ) continue;
2879     if ( !eos._hyp.ToSmooth() ) continue;
2880
2881     const TopoDS_Edge& E = TopoDS::Edge( edgesByGeom[iS]._shape );
2882     if ( SMESH_Algo::isDegenerated( E ) || !edgesOfSmooFaces.Contains( E ))
2883       continue;
2884
2885     double tgtThick = eos._hyp.GetTotalThickness();
2886     for ( TopoDS_Iterator vIt( E ); vIt.More() && !eos._toSmooth; vIt.Next() )
2887     {
2888       TGeomID iV = getMeshDS()->ShapeToIndex( vIt.Value() );
2889       vector<_LayerEdge*>& eV = edgesByGeom[ iV ]._edges;
2890       if ( eV.empty() ) continue;
2891       gp_Vec  eDir    = getEdgeDir( E, TopoDS::Vertex( vIt.Value() ));
2892       double angle    = eDir.Angle( eV[0]->_normal );
2893       double cosin    = Cos( angle );
2894       double cosinAbs = Abs( cosin );
2895       if ( cosinAbs > theMinSmoothCosin )
2896       {
2897         // always smooth analytic EDGEs
2898         Handle(Geom_Curve) curve = _Smoother1D::CurveForSmooth( E, eos, helper );
2899         eos._toSmooth = ! curve.IsNull();
2900
2901         // compare tgtThick with the length of an end segment
2902         SMDS_ElemIteratorPtr eIt = eV[0]->_nodes[0]->GetInverseElementIterator(SMDSAbs_Edge);
2903         while ( eIt->more() && !eos._toSmooth )
2904         {
2905           const SMDS_MeshElement* endSeg = eIt->next();
2906           if ( endSeg->getshapeId() == (int) iS )
2907           {
2908             double segLen =
2909               SMESH_TNodeXYZ( endSeg->GetNode(0) ).Distance( endSeg->GetNode(1 ));
2910             eos._toSmooth = needSmoothing( cosinAbs, tgtThick, segLen );
2911           }
2912         }
2913         if ( eos._toSmooth )
2914         {
2915           eos._edgeSmoother = new _Smoother1D( curve, eos );
2916
2917           for ( size_t i = 0; i < eos._edges.size(); ++i )
2918             eos._edges[i]->Set( _LayerEdge::TO_SMOOTH );
2919         }
2920       }
2921     }
2922     data._nbShapesToSmooth += eos._toSmooth;
2923
2924   } // check EDGEs
2925
2926   // Reset _cosin if no smooth is allowed by the user
2927   for ( size_t iS = 0; iS < edgesByGeom.size(); ++iS )
2928   {
2929     _EdgesOnShape& eos = edgesByGeom[iS];
2930     if ( eos._edges.empty() ) continue;
2931
2932     if ( !eos._hyp.ToSmooth() )
2933       for ( size_t i = 0; i < eos._edges.size(); ++i )
2934         eos._edges[i]->SetCosin( 0 );
2935   }
2936
2937
2938   // Fill _eosC1 to make that C1 FACEs and EGDEs between them to be smoothed as a whole
2939
2940   TopTools_MapOfShape c1VV;
2941
2942   for ( size_t iS = 0; iS < edgesByGeom.size(); ++iS ) // check FACEs
2943   {
2944     _EdgesOnShape& eos = edgesByGeom[iS];
2945     if ( eos._edges.empty() ||
2946          eos.ShapeType() != TopAbs_FACE ||
2947          !eos._toSmooth )
2948       continue;
2949
2950     // check EDGEs of a FACE
2951     TopTools_MapOfShape checkedEE, allVV;
2952     list< SMESH_subMesh* > smQueue( 1, eos._subMesh ); // sm of FACEs
2953     while ( !smQueue.empty() )
2954     {
2955       SMESH_subMesh* sm = smQueue.front();
2956       smQueue.pop_front();
2957       SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(/*includeSelf=*/false);
2958       while ( smIt->more() )
2959       {
2960         sm = smIt->next();
2961         if ( sm->GetSubShape().ShapeType() == TopAbs_VERTEX )
2962           allVV.Add( sm->GetSubShape() );
2963         if ( sm->GetSubShape().ShapeType() != TopAbs_EDGE ||
2964              !checkedEE.Add( sm->GetSubShape() ))
2965           continue;
2966
2967         _EdgesOnShape*      eoe = data.GetShapeEdges( sm->GetId() );
2968         vector<_LayerEdge*>& eE = eoe->_edges;
2969         if ( eE.empty() || !eoe->_sWOL.IsNull() )
2970           continue;
2971
2972         bool isC1 = true; // check continuity along an EDGE
2973         for ( size_t i = 0; i < eE.size() && isC1; ++i )
2974           isC1 = ( Abs( eE[i]->_cosin ) < theMinSmoothCosin );
2975         if ( !isC1 )
2976           continue;
2977
2978         // check that mesh faces are C1 as well
2979         {
2980           gp_XYZ norm1, norm2;
2981           const SMDS_MeshNode*   n = eE[ eE.size() / 2 ]->_nodes[0];
2982           SMDS_ElemIteratorPtr fIt = n->GetInverseElementIterator(SMDSAbs_Face);
2983           if ( !SMESH_MeshAlgos::FaceNormal( fIt->next(), norm1, /*normalized=*/true ))
2984             continue;
2985           while ( fIt->more() && isC1 )
2986             isC1 = ( SMESH_MeshAlgos::FaceNormal( fIt->next(), norm2, /*normalized=*/true ) &&
2987                      Abs( norm1 * norm2 ) >= ( 1. - theMinSmoothCosin ));
2988           if ( !isC1 )
2989             continue;
2990         }
2991
2992         // add the EDGE and an adjacent FACE to _eosC1
2993         PShapeIteratorPtr fIt = helper.GetAncestors( sm->GetSubShape(), *_mesh, TopAbs_FACE );
2994         while ( const TopoDS_Shape* face = fIt->next() )
2995         {
2996           _EdgesOnShape* eof = data.GetShapeEdges( *face );
2997           if ( !eof ) continue; // other solid
2998           if ( !eos.HasC1( eoe ))
2999           {
3000             eos._eosC1.push_back( eoe );
3001             eoe->_toSmooth = false;
3002             data.PrepareEdgesToSmoothOnFace( eoe, /*substituteSrcNodes=*/false );
3003           }
3004           if ( eos._shapeID != eof->_shapeID && !eos.HasC1( eof )) 
3005           {
3006             eos._eosC1.push_back( eof );
3007             eof->_toSmooth = false;
3008             data.PrepareEdgesToSmoothOnFace( eof, /*substituteSrcNodes=*/false );
3009             smQueue.push_back( eof->_subMesh );
3010           }
3011         }
3012       }
3013     }
3014     if ( eos._eosC1.empty() )
3015       continue;
3016
3017     // check VERTEXes of C1 FACEs
3018     TopTools_MapIteratorOfMapOfShape vIt( allVV );
3019     for ( ; vIt.More(); vIt.Next() )
3020     {
3021       _EdgesOnShape* eov = data.GetShapeEdges( vIt.Key() );
3022       if ( !eov || eov->_edges.empty() || !eov->_sWOL.IsNull() )
3023         continue;
3024
3025       bool isC1 = true; // check if all adjacent FACEs are in eos._eosC1
3026       PShapeIteratorPtr fIt = helper.GetAncestors( vIt.Key(), *_mesh, TopAbs_FACE );
3027       while ( const TopoDS_Shape* face = fIt->next() )
3028       {
3029         _EdgesOnShape* eof = data.GetShapeEdges( *face );
3030         if ( !eof ) continue; // other solid
3031         isC1 = ( face->IsSame( eos._shape ) || eos.HasC1( eof ));
3032         if ( !isC1 )
3033           break;
3034       }
3035       if ( isC1 )
3036       {
3037         eos._eosC1.push_back( eov );
3038         data.PrepareEdgesToSmoothOnFace( eov, /*substituteSrcNodes=*/false );
3039         c1VV.Add( eov->_shape );
3040       }
3041     }
3042
3043   } // fill _eosC1 of FACEs
3044
3045
3046   // Find C1 EDGEs
3047
3048   vector< pair< _EdgesOnShape*, gp_XYZ > > dirOfEdges;
3049
3050   for ( size_t iS = 0; iS < edgesByGeom.size(); ++iS ) // check VERTEXes
3051   {
3052     _EdgesOnShape& eov = edgesByGeom[iS];
3053     if ( eov._edges.empty() ||
3054          eov.ShapeType() != TopAbs_VERTEX ||
3055          c1VV.Contains( eov._shape ))
3056       continue;
3057     const TopoDS_Vertex& V = TopoDS::Vertex( eov._shape );
3058
3059     // get directions of surrounding EDGEs
3060     dirOfEdges.clear();
3061     PShapeIteratorPtr fIt = helper.GetAncestors( eov._shape, *_mesh, TopAbs_EDGE );
3062     while ( const TopoDS_Shape* e = fIt->next() )
3063     {
3064       _EdgesOnShape* eoe = data.GetShapeEdges( *e );
3065       if ( !eoe ) continue; // other solid
3066       gp_XYZ eDir = getEdgeDir( TopoDS::Edge( *e ), V );
3067       if ( !Precision::IsInfinite( eDir.X() ))
3068         dirOfEdges.push_back( make_pair( eoe, eDir.Normalized() ));
3069     }
3070
3071     // find EDGEs with C1 directions
3072     for ( size_t i = 0; i < dirOfEdges.size(); ++i )
3073       for ( size_t j = i+1; j < dirOfEdges.size(); ++j )
3074         if ( dirOfEdges[i].first && dirOfEdges[j].first )
3075         {
3076           double dot = dirOfEdges[i].second * dirOfEdges[j].second;
3077           bool isC1 = ( dot < - ( 1. - theMinSmoothCosin ));
3078           if ( isC1 )
3079           {
3080             double maxEdgeLen = 3 * Min( eov._edges[0]->_maxLen, eov._hyp.GetTotalThickness() );
3081             double eLen1 = SMESH_Algo::EdgeLength( TopoDS::Edge( dirOfEdges[i].first->_shape ));
3082             double eLen2 = SMESH_Algo::EdgeLength( TopoDS::Edge( dirOfEdges[j].first->_shape ));
3083             if ( eLen1 < maxEdgeLen ) eov._eosC1.push_back( dirOfEdges[i].first );
3084             if ( eLen2 < maxEdgeLen ) eov._eosC1.push_back( dirOfEdges[j].first );
3085             dirOfEdges[i].first = 0;
3086             dirOfEdges[j].first = 0;
3087           }
3088         }
3089   } // fill _eosC1 of VERTEXes
3090
3091
3092
3093   return ok;
3094 }
3095
3096 //================================================================================
3097 /*!
3098  * \brief initialize data of _EdgesOnShape
3099  */
3100 //================================================================================
3101
3102 void _ViscousBuilder::setShapeData( _EdgesOnShape& eos,
3103                                     SMESH_subMesh* sm,
3104                                     _SolidData&    data )
3105 {
3106   if ( !eos._shape.IsNull() ||
3107        sm->GetSubShape().ShapeType() == TopAbs_WIRE )
3108     return;
3109
3110   SMESH_MesherHelper helper( *_mesh );
3111
3112   eos._subMesh = sm;
3113   eos._shapeID = sm->GetId();
3114   eos._shape   = sm->GetSubShape();
3115   if ( eos.ShapeType() == TopAbs_FACE )
3116     eos._shape.Orientation( helper.GetSubShapeOri( data._solid, eos._shape ));
3117   eos._toSmooth = false;
3118   eos._data = &data;
3119
3120   // set _SWOL
3121   map< TGeomID, TopoDS_Shape >::const_iterator s2s =
3122     data._shrinkShape2Shape.find( eos._shapeID );
3123   if ( s2s != data._shrinkShape2Shape.end() )
3124     eos._sWOL = s2s->second;
3125
3126   eos._isRegularSWOL = true;
3127   if ( eos.SWOLType() == TopAbs_FACE )
3128   {
3129     const TopoDS_Face& F = TopoDS::Face( eos._sWOL );
3130     Handle(ShapeAnalysis_Surface) surface = helper.GetSurface( F );
3131     eos._isRegularSWOL = ( ! surface->HasSingularities( 1e-7 ));
3132   }
3133
3134   // set _hyp
3135   if ( data._hyps.size() == 1 )
3136   {
3137     eos._hyp = data._hyps.back();
3138   }
3139   else
3140   {
3141     // compute average StdMeshers_ViscousLayers parameters
3142     map< TGeomID, const StdMeshers_ViscousLayers* >::iterator f2hyp;
3143     if ( eos.ShapeType() == TopAbs_FACE )
3144     {
3145       if (( f2hyp = data._face2hyp.find( eos._shapeID )) != data._face2hyp.end() )
3146         eos._hyp = f2hyp->second;
3147     }
3148     else
3149     {
3150       PShapeIteratorPtr fIt = helper.GetAncestors( eos._shape, *_mesh, TopAbs_FACE );
3151       while ( const TopoDS_Shape* face = fIt->next() )
3152       {
3153         TGeomID faceID = getMeshDS()->ShapeToIndex( *face );
3154         if (( f2hyp = data._face2hyp.find( faceID )) != data._face2hyp.end() )
3155           eos._hyp.Add( f2hyp->second );
3156       }
3157     }
3158   }
3159
3160   // set _faceNormals
3161   if ( ! eos._hyp.UseSurfaceNormal() )
3162   {
3163     if ( eos.ShapeType() == TopAbs_FACE ) // get normals to elements on a FACE
3164     {
3165       SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
3166       eos._faceNormals.resize( smDS->NbElements() );
3167
3168       SMDS_ElemIteratorPtr eIt = smDS->GetElements();
3169       for ( int iF = 0; eIt->more(); ++iF )
3170       {
3171         const SMDS_MeshElement* face = eIt->next();
3172         if ( !SMESH_MeshAlgos::FaceNormal( face, eos._faceNormals[iF], /*normalized=*/true ))
3173           eos._faceNormals[iF].SetCoord( 0,0,0 );
3174       }
3175
3176       if ( !helper.IsReversedSubMesh( TopoDS::Face( eos._shape )))
3177         for ( size_t iF = 0; iF < eos._faceNormals.size(); ++iF )
3178           eos._faceNormals[iF].Reverse();
3179     }
3180     else // find EOS of adjacent FACEs
3181     {
3182       PShapeIteratorPtr fIt = helper.GetAncestors( eos._shape, *_mesh, TopAbs_FACE );
3183       while ( const TopoDS_Shape* face = fIt->next() )
3184       {
3185         TGeomID faceID = getMeshDS()->ShapeToIndex( *face );
3186         eos._faceEOS.push_back( & data._edgesOnShape[ faceID ]);
3187         if ( eos._faceEOS.back()->_shape.IsNull() )
3188           // avoid using uninitialised _shapeID in GetNormal()
3189           eos._faceEOS.back()->_shapeID = faceID;
3190       }
3191     }
3192   }
3193 }
3194
3195 //================================================================================
3196 /*!
3197  * \brief Returns normal of a face
3198  */
3199 //================================================================================
3200
3201 bool _EdgesOnShape::GetNormal( const SMDS_MeshElement* face, gp_Vec& norm )
3202 {
3203   bool ok = false;
3204   const _EdgesOnShape* eos = 0;
3205
3206   if ( face->getshapeId() == _shapeID )
3207   {
3208     eos = this;
3209   }
3210   else
3211   {
3212     for ( size_t iF = 0; iF < _faceEOS.size() && !eos; ++iF )
3213       if ( face->getshapeId() == _faceEOS[ iF ]->_shapeID )
3214         eos = _faceEOS[ iF ];
3215   }
3216
3217   if (( eos ) &&
3218       ( ok = ( face->getIdInShape() < (int) eos->_faceNormals.size() )))
3219   {
3220     norm = eos->_faceNormals[ face->getIdInShape() ];
3221   }
3222   else if ( !eos )
3223   {
3224     debugMsg( "_EdgesOnShape::Normal() failed for face "<<face->GetID()
3225               << " on _shape #" << _shapeID );
3226   }
3227   return ok;
3228 }
3229
3230
3231 //================================================================================
3232 /*!
3233  * \brief Set data of _LayerEdge needed for smoothing
3234  */
3235 //================================================================================
3236
3237 bool _ViscousBuilder::setEdgeData(_LayerEdge&         edge,
3238                                   _EdgesOnShape&      eos,
3239                                   SMESH_MesherHelper& helper,
3240                                   _SolidData&         data)
3241 {
3242   const SMDS_MeshNode* node = edge._nodes[0]; // source node
3243
3244   edge._len       = 0;
3245   edge._maxLen    = Precision::Infinite();
3246   edge._minAngle  = 0;
3247   edge._2neibors  = 0;
3248   edge._curvature = 0;
3249   edge._flags     = 0;
3250
3251   // --------------------------
3252   // Compute _normal and _cosin
3253   // --------------------------
3254
3255   edge._cosin     = 0;
3256   edge._lenFactor = 1.;
3257   edge._normal.SetCoord(0,0,0);
3258   _Simplex::GetSimplices( node, edge._simplices, data._ignoreFaceIds, &data );
3259
3260   int totalNbFaces = 0;
3261   TopoDS_Face F;
3262   std::pair< TopoDS_Face, gp_XYZ > face2Norm[20];
3263   gp_Vec geomNorm;
3264   bool normOK = true;
3265
3266   const bool onShrinkShape = !eos._sWOL.IsNull();
3267   const bool useGeometry   = (( eos._hyp.UseSurfaceNormal() ) ||
3268                               ( eos.ShapeType() != TopAbs_FACE /*&& !onShrinkShape*/ ));
3269
3270   // get geom FACEs the node lies on
3271   //if ( useGeometry )
3272   {
3273     set<TGeomID> faceIds;
3274     if  ( eos.ShapeType() == TopAbs_FACE )
3275     {
3276       faceIds.insert( eos._shapeID );
3277     }
3278     else
3279     {
3280       SMDS_ElemIteratorPtr fIt = node->GetInverseElementIterator(SMDSAbs_Face);
3281       while ( fIt->more() )
3282         faceIds.insert( fIt->next()->getshapeId() );
3283     }
3284     set<TGeomID>::iterator id = faceIds.begin();
3285     for ( ; id != faceIds.end(); ++id )
3286     {
3287       const TopoDS_Shape& s = getMeshDS()->IndexToShape( *id );
3288       if ( s.IsNull() || s.ShapeType() != TopAbs_FACE || data._ignoreFaceIds.count( *id ))
3289         continue;
3290       F = TopoDS::Face( s );
3291       face2Norm[ totalNbFaces ].first = F;
3292       totalNbFaces++;
3293     }