1 // Copyright (C) 2007-2016 CEA/DEN, EDF R&D, OPEN CASCADE
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.
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.
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
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
20 // File : StdMeshers_ViscousLayers.cxx
21 // Created : Wed Dec 1 15:15:34 2010
22 // Author : Edward AGAPOV (eap)
24 #include "StdMeshers_ViscousLayers.hxx"
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"
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>
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>
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>
80 #include <TopoDS_Edge.hxx>
81 #include <TopoDS_Face.hxx>
82 #include <TopoDS_Vertex.hxx>
84 #include <gp_Cone.hxx>
85 #include <gp_Sphere.hxx>
97 //#define __NOT_INVALIDATE_BAD_SMOOTH
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
106 //================================================================================
111 enum UIndex { U_TGT = 1, U_SRC, LEN_TGT };
113 const double theMinSmoothCosin = 0.1;
114 const double theSmoothThickToElemSizeRatio = 0.3;
115 const double theMinSmoothTriaAngle = 30;
116 const double theMinSmoothQuadAngle = 45;
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;
122 bool needSmoothing( double cosin, double tgtThick, double elemSize )
124 return cosin * tgtThick > theSmoothThickToElemSizeRatio * elemSize;
126 double getSmoothingThickness( double cosin, double elemSize )
128 return theSmoothThickToElemSizeRatio * elemSize / cosin;
132 * \brief SMESH_ProxyMesh computed by _ViscousBuilder for a SOLID.
133 * It is stored in a SMESH_subMesh of the SOLID as SMESH_subMeshEventListenerData
135 struct _MeshOfSolid : public SMESH_ProxyMesh,
136 public SMESH_subMeshEventListenerData
138 bool _n2nMapComputed;
139 SMESH_ComputeErrorPtr _warning;
141 _MeshOfSolid( SMESH_Mesh* mesh)
142 :SMESH_subMeshEventListenerData( /*isDeletable=*/true),_n2nMapComputed(false)
144 SMESH_ProxyMesh::setMesh( *mesh );
147 // returns submesh for a geom face
148 SMESH_ProxyMesh::SubMesh* getFaceSubM(const TopoDS_Face& F, bool create=false)
150 TGeomID i = SMESH_ProxyMesh::shapeIndex(F);
151 return create ? SMESH_ProxyMesh::getProxySubMesh(i) : findProxySubMesh(i);
153 void setNode2Node(const SMDS_MeshNode* srcNode,
154 const SMDS_MeshNode* proxyNode,
155 const SMESH_ProxyMesh::SubMesh* subMesh)
157 SMESH_ProxyMesh::setNode2Node( srcNode,proxyNode,subMesh);
160 //--------------------------------------------------------------------------------
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
165 class _ShrinkShapeListener : SMESH_subMeshEventListener
167 _ShrinkShapeListener()
168 : SMESH_subMeshEventListener(/*isDeletable=*/false,
169 "StdMeshers_ViscousLayers::_ShrinkShapeListener") {}
171 static SMESH_subMeshEventListener* Get() { static _ShrinkShapeListener l; return &l; }
172 virtual void ProcessEvent(const int event,
174 SMESH_subMesh* solidSM,
175 SMESH_subMeshEventListenerData* data,
176 const SMESH_Hypothesis* hyp)
178 if ( SMESH_subMesh::COMPUTE_EVENT == eventType && solidSM->IsEmpty() && data )
180 SMESH_subMeshEventListener::ProcessEvent(event,eventType,solidSM,data,hyp);
184 //--------------------------------------------------------------------------------
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
190 class _ViscousListener : SMESH_subMeshEventListener
193 SMESH_subMeshEventListener(/*isDeletable=*/false,
194 "StdMeshers_ViscousLayers::_ViscousListener") {}
195 static SMESH_subMeshEventListener* Get() { static _ViscousListener l; return &l; }
197 virtual void ProcessEvent(const int event,
199 SMESH_subMesh* subMesh,
200 SMESH_subMeshEventListenerData* data,
201 const SMESH_Hypothesis* hyp)
203 if (( SMESH_subMesh::COMPUTE_EVENT == eventType ) &&
204 ( SMESH_subMesh::CHECK_COMPUTE_STATE != event &&
205 SMESH_subMesh::SUBMESH_COMPUTED != event ))
207 // delete SMESH_ProxyMesh containing temporary faces
208 subMesh->DeleteEventListener( this );
211 // Finds or creates proxy mesh of the solid
212 static _MeshOfSolid* GetSolidMesh(SMESH_Mesh* mesh,
213 const TopoDS_Shape& solid,
216 if ( !mesh ) return 0;
217 SMESH_subMesh* sm = mesh->GetSubMesh(solid);
218 _MeshOfSolid* data = (_MeshOfSolid*) sm->GetEventListenerData( Get() );
219 if ( !data && toCreate )
221 data = new _MeshOfSolid(mesh);
222 data->mySubMeshes.push_back( sm ); // to find SOLID by _MeshOfSolid
223 sm->SetEventListener( Get(), data, sm );
227 // Removes proxy mesh of the solid
228 static void RemoveSolidMesh(SMESH_Mesh* mesh, const TopoDS_Shape& solid)
230 mesh->GetSubMesh(solid)->DeleteEventListener( _ViscousListener::Get() );
234 //================================================================================
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
241 //================================================================================
243 void ToClearSubWithMain( SMESH_subMesh* sub, const TopoDS_Shape& main)
245 SMESH_subMesh* mainSM = sub->GetFather()->GetSubMesh( main );
246 SMESH_subMeshEventListenerData* data =
247 mainSM->GetEventListenerData( _ShrinkShapeListener::Get());
250 if ( find( data->mySubMeshes.begin(), data->mySubMeshes.end(), sub ) ==
251 data->mySubMeshes.end())
252 data->mySubMeshes.push_back( sub );
256 data = SMESH_subMeshEventListenerData::MakeData( /*dependent=*/sub );
257 sub->SetEventListener( _ShrinkShapeListener::Get(), data, /*whereToListenTo=*/mainSM );
261 //--------------------------------------------------------------------------------
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.
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
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]);
290 bool IsForward(const SMDS_MeshNode* nSrc, const gp_XYZ& pTgt, double& vol) const
292 SMESH_TNodeXYZ pSrc( nSrc );
293 return IsForward( &pSrc, &pTgt, vol );
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
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 );
305 return d*refSign > 1e-100;
307 bool IsMinAngleOK( const gp_XYZ& pTgt, double& minAngle ) const
309 SMESH_TNodeXYZ pPrev( _nPrev ), pNext( _nNext );
310 if ( !_nOpp ) // triangle
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();
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;
322 minAngle = ( pn * -tp ) * ( pn * -tp ) / pn2 / tp2;
324 static double theMaxCos2 = ( Cos( theMinSmoothTriaAngle * M_PI / 180. ) *
325 Cos( theMinSmoothTriaAngle * M_PI / 180. ));
326 return minAngle < theMaxCos2;
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 )));
341 static double theMaxCos2 = ( Cos( theMinSmoothQuadAngle * M_PI / 180. ) *
342 Cos( theMinSmoothQuadAngle * M_PI / 180. ));
343 return minAngle < theMaxCos2;
346 bool IsNeighbour(const _Simplex& other) const
348 return _nPrev == other._nNext || _nNext == other._nPrev;
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);
358 //--------------------------------------------------------------------------------
360 * Structure used to take into account surface curvature while smoothing
365 double _k; // factor to correct node smoothed position
366 double _h2lenRatio; // avgNormProj / (2*avgDist)
367 gp_Pnt2d _uv; // UV used in putOnOffsetSurface()
369 static _Curvature* New( double avgNormProj, double avgDist )
372 if ( fabs( avgNormProj / avgDist ) > 1./200 )
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 );
381 c->_uv.SetCoord( 0., 0. );
385 double lenDelta(double len) const { return _k * ( _r + len ); }
386 double lenDeltaByDist(double dist) const { return dist * _h2lenRatio; }
388 //--------------------------------------------------------------------------------
392 struct _EdgesOnShape;
394 typedef map< const SMDS_MeshNode*, _LayerEdge*, TIDCompare > TNode2Edge;
396 //--------------------------------------------------------------------------------
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())
403 typedef gp_XYZ (_LayerEdge::*PSmooFun)();
405 vector< const SMDS_MeshNode*> _nodes;
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
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;
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
440 bool Is ( EFlags f ) const { return _flags & f; }
441 void Set ( EFlags f ) { _flags |= f; }
442 void Unset( EFlags f ) { _flags &= ~f; }
444 void SetNewLength( double len, _EdgesOnShape& eos, SMESH_MesherHelper& helper );
445 bool SetNewLength2d( Handle(Geom_Surface)& surface,
446 const TopoDS_Face& F,
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,
468 vector< _LayerEdge* > & badSmooEdges);
469 bool FindIntersection( SMESH_ElementSearcher& searcher,
471 const double& epsilon,
473 const SMDS_MeshElement** face = 0);
474 bool SegTriaInter( const gp_Ax1& lastSegment,
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,
485 const double& epsilon) const
486 { return SegTriaInter( lastSegment,
487 SMESH_TNodeXYZ( n0 ), SMESH_TNodeXYZ( n1 ), SMESH_TNodeXYZ( n2 ),
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
503 double GetSmooLen() { return _cosin; } // for _LayerEdge's on FACE _cosin is not used
505 gp_XYZ smoothLaplacian();
506 gp_XYZ smoothAngular();
507 gp_XYZ smoothLengthWeighted();
508 gp_XYZ smoothCentroidal();
509 gp_XYZ smoothNefPolygon();
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;
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",
530 bool operator () (const _LayerEdge* e1, const _LayerEdge* e2) const
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 );
536 //--------------------------------------------------------------------------------
538 * A 2D half plane used by _LayerEdge::smoothNefPolygon()
542 gp_XY _pos, _dir, _inNorm;
543 bool IsOut( const gp_XY p, const double tol ) const
545 return _inNorm * ( p - _pos ) < -tol;
547 bool FindIntersection( const _halfPlane& hp, gp_XY & intPnt )
549 //const double eps = 1e-10;
550 double D = _dir.Crossed( hp._dir );
551 if ( fabs(D) < std::numeric_limits<double>::min())
553 gp_XY vec21 = _pos - hp._pos;
554 double u = hp._dir.Crossed( vec21 ) / D;
555 intPnt = _pos + _dir * u;
559 //--------------------------------------------------------------------------------
561 * Structure used to smooth a _LayerEdge based on an EDGE.
565 double _wgt [2]; // weights of _nodes
566 _LayerEdge* _edges[2];
568 // normal to plane passing through _LayerEdge._normal and tangent of EDGE
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;
575 const SMDS_MeshNode* srcNode(bool is2nd) {
576 return _edges[is2nd] ? _edges[is2nd]->_nodes[0] : 0;
579 std::swap( _wgt [0], _wgt [1] );
580 std::swap( _edges[0], _edges[1] );
582 void set( _LayerEdge* e1, _LayerEdge* e2, double w1, double w2 ) {
583 _edges[0] = e1; _edges[1] = e2; _wgt[0] = w1; _wgt[1] = w2;
585 bool include( const _LayerEdge* e ) {
586 return ( _edges[0] == e || _edges[1] == e );
591 //--------------------------------------------------------------------------------
593 * \brief Layers parameters got by averaging several hypotheses
597 AverageHyp( const StdMeshers_ViscousLayers* hyp = 0 )
598 :_nbLayers(0), _nbHyps(0), _method(0), _thickness(0), _stretchFactor(0)
602 void Add( const StdMeshers_ViscousLayers* hyp )
607 _nbLayers = hyp->GetNumberLayers();
608 //_thickness += hyp->GetTotalThickness();
609 _thickness = Max( _thickness, hyp->GetTotalThickness() );
610 _stretchFactor += hyp->GetStretchFactor();
611 _method = hyp->GetMethod();
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; }
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; }
627 int _nbLayers, _nbHyps, _method;
628 double _thickness, _stretchFactor;
631 //--------------------------------------------------------------------------------
633 * \brief _LayerEdge's on a shape and other shape data
637 vector< _LayerEdge* > _edges;
641 SMESH_subMesh * _subMesh;
642 // face or edge w/o layer along or near which _edges are inflated
644 bool _isRegularSWOL; // w/o singularities
645 // averaged StdMeshers_ViscousLayers parameters
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
652 vector< gp_XYZ > _faceNormals; // if _shape is FACE
653 vector< _EdgesOnShape* > _faceEOS; // to get _faceNormals of adjacent FACEs
655 Handle(ShapeAnalysis_Surface) _offsetSurf;
656 _LayerEdge* _edgeForOffset;
658 _SolidData* _data; // parent SOLID
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; }
669 _EdgesOnShape(): _shapeID(-1), _subMesh(0), _toSmooth(false), _edgeSmoother(0) {}
672 //--------------------------------------------------------------------------------
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
682 // edges whose _simplices are used to detect prism distortion
683 vector< _LayerEdge* > _simplexTestEdges;
685 // map a sub-shape to _SolidData::_edgesOnShape
686 map< TGeomID, _EdgesOnShape* > _subIdToEOS;
690 bool GetCenterOfCurvature( _LayerEdge* ledge,
691 BRepLProp_SLProps& surfProp,
692 SMESH_MesherHelper& helper,
693 gp_Pnt & center ) const;
694 bool CheckPrisms() const;
697 //--------------------------------------------------------------------------------
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()
703 struct _CollisionEdges
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(); }
711 //--------------------------------------------------------------------------------
713 * \brief Data of a SOLID
717 typedef const StdMeshers_ViscousLayers* THyp;
719 TGeomID _index; // SOLID id
720 _MeshOfSolid* _proxyMesh;
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
727 double _stepSize, _stepSizeCoeff, _geomSize;
728 const SMDS_MeshNode* _stepSizeNodes[2];
730 TNode2Edge _n2eMap; // nodes and _LayerEdge's based on them
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;
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;
743 // Convex FACEs whose radius of curvature is less than the thickness of layers
744 map< TGeomID, _ConvexFace > _convexFaces;
746 // shapes (EDGEs and VERTEXes) srink from which is forbidden due to collisions with
747 // the adjacent SOLID
748 set< TGeomID > _noShrinkShapes;
750 int _nbShapesToSmooth;
752 //map< TGeomID,Handle(Geom_Curve)> _edge2curve;
754 vector< _CollisionEdges > _collisionEdges;
755 set< TGeomID > _concaveFaces;
757 double _maxThickness; // of all _hyps
758 double _minThickness; // of all _hyps
760 double _epsilon; // precision for SegTriaInter()
762 SMESH_MesherHelper* _helper;
764 _SolidData(const TopoDS_Shape& s=TopoDS_Shape(),
766 :_solid(s), _proxyMesh(m), _helper(0) {}
769 void SortOnEdge( const TopoDS_Edge& E, vector< _LayerEdge* >& edges);
770 void Sort2NeiborsOnEdge( vector< _LayerEdge* >& edges );
772 _ConvexFace* GetConvexFace( const TGeomID faceID ) {
773 map< TGeomID, _ConvexFace >::iterator id2face = _convexFaces.find( faceID );
774 return id2face == _convexFaces.end() ? 0 : & id2face->second;
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() ); }
781 SMESH_MesherHelper& GetHelper() const { return *_helper; }
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 );
788 void AddShapesToSmooth( const set< _EdgesOnShape* >& shape,
789 const set< _EdgesOnShape* >* edgesNoAnaSmooth=0 );
791 void PrepareEdgesToSmoothOnFace( _EdgesOnShape* eof, bool substituteSrcNodes );
793 //--------------------------------------------------------------------------------
795 * \brief Offset plane used in getNormalByOffset()
801 int _faceIndexNext[2];
802 gp_Lin _lines[2]; // line of intersection with neighbor _OffsetPlane's
805 _isLineOK[0] = _isLineOK[1] = false; _faceIndexNext[0] = _faceIndexNext[1] = -1;
807 void ComputeIntersectionLine( _OffsetPlane& pln );
808 gp_XYZ GetCommonPoint(bool& isFound) const;
809 int NbLines() const { return _isLineOK[0] + _isLineOK[1]; }
811 //--------------------------------------------------------------------------------
813 * \brief Container of centers of curvature at nodes on an EDGE bounding _ConvexFace
815 struct _CentralCurveOnEdge
818 vector< gp_Pnt > _curvaCenters;
819 vector< _LayerEdge* > _ledges;
820 vector< gp_XYZ > _normals; // new normal for each of _ledges
821 vector< double > _segLength2;
824 TopoDS_Face _adjFace;
825 bool _adjFaceToSmooth;
827 void Append( const gp_Pnt& center, _LayerEdge* ledge )
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 );
835 bool FindNewNormal( const gp_Pnt& center, gp_XYZ& newNormal );
836 void SetShapes( const TopoDS_Edge& edge,
837 const _ConvexFace& convFace,
839 SMESH_MesherHelper& helper);
841 //--------------------------------------------------------------------------------
843 * \brief Data of node on a shrinked FACE
847 const SMDS_MeshNode* _node;
848 vector<_Simplex> _simplices; // for quality check
850 enum SmoothType { LAPLACIAN, CENTROIDAL, ANGULAR, TFI };
852 bool Smooth(int& badNb,
853 Handle(Geom_Surface)& surface,
854 SMESH_MesherHelper& helper,
855 const double refSign,
859 gp_XY computeAngularPos(vector<gp_XY>& uv,
860 const gp_XY& uvToFix,
861 const double refSign );
863 //--------------------------------------------------------------------------------
865 * \brief Builder of viscous layers
867 class _ViscousBuilder
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 );
878 // restore event listeners used to clear an inferior dim sub-mesh modified by viscous layers
879 void RestoreListeners();
881 // computes SMESH_ProxyMesh::SubMesh::_n2n;
882 bool MakeN2NMap( _MeshOfSolid* pm );
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,
900 bool shiftInside=false);
901 bool getFaceNormalAtSingularity(const gp_XY& uv,
902 const TopoDS_Face& face,
903 SMESH_MesherHelper& helper,
905 gp_XYZ getWeigthedNormal( const _LayerEdge* edge );
906 gp_XYZ getNormalByOffset( _LayerEdge* edge,
907 std::pair< TopoDS_Face, gp_XYZ > fId2Normal[],
909 bool findNeiborsOnEdge(const _LayerEdge* edge,
910 const SMDS_MeshNode*& n1,
911 const SMDS_MeshNode*& n2,
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,
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,
937 void updateNormalsOfC1Vertices( _SolidData& data );
938 bool updateNormalsOfSmoothed( _SolidData& data,
939 SMESH_MesherHelper& helper,
941 const double stepSize );
942 bool isNewNormalOk( _SolidData& data,
944 const gp_XYZ& newNormal);
945 bool refine(_SolidData& data);
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,
955 set<const SMDS_MeshNode*> * involvedNodes=NULL);
956 bool addBoundaryElements();
958 bool error( const string& text, int solidID=-1 );
959 SMESHDS_Mesh* getMeshDS() const { return _mesh->GetMeshDS(); }
962 void makeGroupOfLE();
965 SMESH_ComputeErrorPtr _error;
967 vector< _SolidData > _sdVec;
970 //--------------------------------------------------------------------------------
972 * \brief Shrinker of nodes on the EDGE
976 TopoDS_Edge _geomEdge;
977 vector<double> _initU;
978 vector<double> _normPar;
979 vector<const SMDS_MeshNode*> _nodes;
980 const _LayerEdge* _edges[2];
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; }
993 //--------------------------------------------------------------------------------
995 * \brief Smoother of _LayerEdge's on EDGE.
999 struct OffPnt // point of the offsetted EDGE
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(); }
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;
1013 static Handle(Geom_Curve) CurveForSmooth( const TopoDS_Edge& E,
1015 SMESH_MesherHelper& helper);
1017 _Smoother1D( Handle(Geom_Curve) curveForSmooth,
1018 _EdgesOnShape& eos )
1019 : _anaCurve( curveForSmooth ), _eos( eos )
1022 bool Perform(_SolidData& data,
1023 Handle(ShapeAnalysis_Surface)& surface,
1024 const TopoDS_Face& F,
1025 SMESH_MesherHelper& helper )
1027 if ( _leParams.empty() || ( !isAnalytic() && _offPoints.empty() ))
1031 return smoothAnalyticEdge( data, surface, F, helper );
1033 return smoothComplexEdge ( data, surface, F, helper );
1035 void prepare(_SolidData& data );
1037 bool smoothAnalyticEdge( _SolidData& data,
1038 Handle(ShapeAnalysis_Surface)& surface,
1039 const TopoDS_Face& F,
1040 SMESH_MesherHelper& helper);
1042 bool smoothComplexEdge( _SolidData& data,
1043 Handle(ShapeAnalysis_Surface)& surface,
1044 const TopoDS_Face& F,
1045 SMESH_MesherHelper& helper);
1047 void setNormalOnV( const bool is2nd,
1048 SMESH_MesherHelper& helper);
1050 _LayerEdge* getLEdgeOnV( bool is2nd )
1052 return _eos._edges[ is2nd ? _eos._edges.size()-1 : 0 ]->_2neibors->_edges[ is2nd ];
1054 bool isAnalytic() const { return !_anaCurve.IsNull(); }
1056 //--------------------------------------------------------------------------------
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
1062 struct _TmpMeshFace : public SMDS_MeshElement
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()));}
1077 //--------------------------------------------------------------------------------
1079 * \brief Class of temporary mesh face storing _LayerEdge it's based on
1081 struct _TmpMeshFaceOnEdge : public _TmpMeshFace
1083 _LayerEdge *_le1, *_le2;
1084 _TmpMeshFaceOnEdge( _LayerEdge* le1, _LayerEdge* le2, int ID ):
1085 _TmpMeshFace( vector<const SMDS_MeshNode*>(4), ID ), _le1(le1), _le2(le2)
1087 _nn[0]=_le1->_nodes[0];
1088 _nn[1]=_le1->_nodes.back();
1089 _nn[2]=_le2->_nodes.back();
1090 _nn[3]=_le2->_nodes[0];
1092 gp_XYZ GetDir() const // return average direction of _LayerEdge's, normal to EDGE
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 );
1106 gp_XYZ GetDir(_LayerEdge* le1, _LayerEdge* le2) // return average direction of _LayerEdge's
1108 _nn[0]=le1->_nodes[0];
1109 _nn[1]=le1->_nodes.back();
1110 _nn[2]=le2->_nodes.back();
1111 _nn[3]=le2->_nodes[0];
1115 //--------------------------------------------------------------------------------
1117 * \brief Retriever of node coordinates either directly or from a surface by node UV.
1118 * \warning Location of a surface is ignored
1120 struct _NodeCoordHelper
1122 SMESH_MesherHelper& _helper;
1123 const TopoDS_Face& _face;
1124 Handle(Geom_Surface) _surface;
1125 gp_XYZ (_NodeCoordHelper::* _fun)(const SMDS_MeshNode* n) const;
1127 _NodeCoordHelper(const TopoDS_Face& F, SMESH_MesherHelper& helper, bool is2D)
1128 : _helper( helper ), _face( F )
1132 TopLoc_Location loc;
1133 _surface = BRep_Tool::Surface( _face, loc );
1135 if ( _surface.IsNull() )
1136 _fun = & _NodeCoordHelper::direct;
1138 _fun = & _NodeCoordHelper::byUV;
1140 gp_XYZ operator()(const SMDS_MeshNode* n) const { return (this->*_fun)( n ); }
1143 gp_XYZ direct(const SMDS_MeshNode* n) const
1145 return SMESH_TNodeXYZ( n );
1147 gp_XYZ byUV (const SMDS_MeshNode* n) const
1149 gp_XY uv = _helper.GetNodeUV( _face, n );
1150 return _surface->Value( uv.X(), uv.Y() ).XYZ();
1154 //================================================================================
1156 * \brief Check angle between vectors
1158 //================================================================================
1160 inline bool isLessAngle( const gp_Vec& v1, const gp_Vec& v2, const double cos )
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 ));
1169 } // namespace VISCOUS_3D
1173 //================================================================================
1174 // StdMeshers_ViscousLayers hypothesis
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 )
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)
1186 if ( faceIds != _shapeIds )
1187 _shapeIds = faceIds, NotifySubMeshesHypothesisModification();
1188 if ( _isToIgnoreShapes != toIgnore )
1189 _isToIgnoreShapes = toIgnore, NotifySubMeshesHypothesisModification();
1190 } // --------------------------------------------------------------------------------
1191 void StdMeshers_ViscousLayers::SetTotalThickness(double thickness)
1193 if ( thickness != _thickness )
1194 _thickness = thickness, NotifySubMeshesHypothesisModification();
1195 } // --------------------------------------------------------------------------------
1196 void StdMeshers_ViscousLayers::SetNumberLayers(int nb)
1198 if ( _nbLayers != nb )
1199 _nbLayers = nb, NotifySubMeshesHypothesisModification();
1200 } // --------------------------------------------------------------------------------
1201 void StdMeshers_ViscousLayers::SetStretchFactor(double factor)
1203 if ( _stretchFactor != factor )
1204 _stretchFactor = factor, NotifySubMeshesHypothesisModification();
1205 } // --------------------------------------------------------------------------------
1206 void StdMeshers_ViscousLayers::SetMethod( ExtrusionMethod method )
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
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();
1222 vector<SMESH_ProxyMesh::Ptr> components;
1223 TopExp_Explorer exp( theShape, TopAbs_SOLID );
1224 for ( ; exp.More(); exp.Next() )
1226 if ( _MeshOfSolid* pm =
1227 _ViscousListener::GetSolidMesh( &theMesh, exp.Current(), /*toCreate=*/false))
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
1235 if ( pm->_warning && !pm->_warning->IsOK() )
1237 SMESH_subMesh* sm = theMesh.GetSubMesh( exp.Current() );
1238 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1239 if ( !smError || smError->IsOK() )
1240 smError = pm->_warning;
1243 _ViscousListener::RemoveSolidMesh ( &theMesh, exp.Current() );
1245 switch ( components.size() )
1249 case 1: return components[0];
1251 default: return SMESH_ProxyMesh::Ptr( new SMESH_ProxyMesh( components ));
1253 return SMESH_ProxyMesh::Ptr();
1254 } // --------------------------------------------------------------------------------
1255 std::ostream & StdMeshers_ViscousLayers::SaveTo(std::ostream & save)
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;
1266 } // --------------------------------------------------------------------------------
1267 std::istream & StdMeshers_ViscousLayers::LoadFrom(std::istream & load)
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;
1279 _isToIgnoreShapes = true; // old behavior
1282 } // --------------------------------------------------------------------------------
1283 bool StdMeshers_ViscousLayers::SetParametersByMesh(const SMESH_Mesh* theMesh,
1284 const TopoDS_Shape& theShape)
1288 } // --------------------------------------------------------------------------------
1289 SMESH_ComputeErrorPtr
1290 StdMeshers_ViscousLayers::CheckHypothesis(SMESH_Mesh& theMesh,
1291 const TopoDS_Shape& theShape,
1292 SMESH_Hypothesis::Hypothesis_Status& theStatus)
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;
1299 theStatus = SMESH_Hypothesis::HYP_OK;
1303 // --------------------------------------------------------------------------------
1304 bool StdMeshers_ViscousLayers::IsShapeWithLayers(int shapeIndex) const
1307 ( std::find( _shapeIds.begin(), _shapeIds.end(), shapeIndex ) != _shapeIds.end() );
1308 return IsToIgnoreShapes() ? !isIn : isIn;
1310 // END StdMeshers_ViscousLayers hypothesis
1311 //================================================================================
1313 namespace VISCOUS_3D
1315 gp_XYZ getEdgeDir( const TopoDS_Edge& E, const TopoDS_Vertex& fromV )
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();
1328 //--------------------------------------------------------------------------------
1329 gp_XYZ getEdgeDir( const TopoDS_Edge& E, const SMDS_MeshNode* atNode,
1330 SMESH_MesherHelper& helper)
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 );
1340 //--------------------------------------------------------------------------------
1341 gp_XYZ getFaceDir( const TopoDS_Face& F, const TopoDS_Vertex& fromV,
1342 const SMDS_MeshNode* node, SMESH_MesherHelper& helper, bool& ok,
1344 //--------------------------------------------------------------------------------
1345 gp_XYZ getFaceDir( const TopoDS_Face& F, const TopoDS_Edge& fromE,
1346 const SMDS_MeshNode* node, SMESH_MesherHelper& helper, bool& ok)
1349 Handle(Geom_Curve) c = BRep_Tool::Curve( fromE, f, l );
1352 TopoDS_Vertex v = helper.IthVertex( 0, fromE );
1353 return getFaceDir( F, v, node, helper, ok );
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 );
1361 double u = helper.GetNodeU( fromE, node, 0, &ok );
1363 TopAbs_Orientation o = helper.GetSubShapeOri( F.Oriented(TopAbs_FORWARD), fromE);
1364 if ( o == TopAbs_REVERSED )
1367 gp_Vec dir = norm ^ du;
1369 if ( node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX &&
1370 helper.IsClosedEdge( fromE ))
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 )
1376 gp_Vec dir2 = norm ^ dv;
1377 dir = dir.Normalized() + dir2.Normalized();
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)
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
1392 TopoDS_Vertex VV[2];
1393 TopExp_Explorer exp( faceFrw, TopAbs_EDGE );
1394 for ( ; exp.More() && nbEdges < 2; exp.Next() )
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();
1403 else if ( VV[0].IsSame( fromV )) {
1404 nbEdges += edges[ 1 ].IsNull();
1409 gp_XYZ dir(0,0,0), edgeDir[2];
1412 // get dirs of edges going fromV
1414 for ( size_t i = 0; i < nbEdges && ok; ++i )
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 );
1421 if ( !ok ) return dir;
1423 // get angle between the 2 edges
1425 double angle = helper.GetAngle( edges[0], edges[1], faceFrw, fromV, &faceNormal );
1426 if ( Abs( angle ) < 5 * M_PI/180 )
1428 dir = ( faceNormal.XYZ() ^ edgeDir[0].Reversed()) + ( faceNormal.XYZ() ^ edgeDir[1] );
1432 dir = edgeDir[0] + edgeDir[1];
1437 double angle = gp_Vec( edgeDir[0] ).Angle( dir );
1438 *cosin = Cos( angle );
1441 else if ( nbEdges == 1 )
1443 dir = getFaceDir( faceFrw, edges[ edges[0].IsNull() ], node, helper, ok );
1444 if ( cosin ) *cosin = 1.;
1454 //================================================================================
1456 * \brief Finds concave VERTEXes of a FACE
1458 //================================================================================
1460 bool getConcaveVertices( const TopoDS_Face& F,
1461 SMESH_MesherHelper& helper,
1462 set< TGeomID >* vertices = 0)
1464 // check angles at VERTEXes
1466 TSideVector wires = StdMeshers_FaceSide::GetFaceWires( F, *helper.GetMesh(), 0, error );
1467 for ( size_t iW = 0; iW < wires.size(); ++iW )
1469 const int nbEdges = wires[iW]->NbEdges();
1470 if ( nbEdges < 2 && SMESH_Algo::isDegenerated( wires[iW]->Edge(0)))
1472 for ( int iE1 = 0; iE1 < nbEdges; ++iE1 )
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. )
1485 vertices->insert( helper.GetMeshDS()->ShapeToIndex( V ));
1489 return vertices ? !vertices->empty() : false;
1492 //================================================================================
1494 * \brief Returns true if a FACE is bound by a concave EDGE
1496 //================================================================================
1498 bool isConcave( const TopoDS_Face& F,
1499 SMESH_MesherHelper& helper,
1500 set< TGeomID >* vertices = 0 )
1502 bool isConcv = false;
1503 // if ( helper.Count( F, TopAbs_WIRE, /*useMap=*/false) > 1 )
1505 gp_Vec2d drv1, drv2;
1507 TopExp_Explorer eExp( F.Oriented( TopAbs_FORWARD ), TopAbs_EDGE );
1508 for ( ; eExp.More(); eExp.Next() )
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 )
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 )
1526 isConvex = ( cross > -1e-9 ); // 0.1 );
1530 //cout << "Concave FACE " << helper.GetMeshDS()->ShapeToIndex( F ) << endl;
1539 // check angles at VERTEXes
1540 if ( getConcaveVertices( F, helper, vertices ))
1546 //================================================================================
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
1554 //================================================================================
1556 bool getDistFromEdge( const SMDS_MeshElement* face,
1557 const SMDS_MeshNode* nodeOnEdge,
1560 faceSize = Precision::Infinite();
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 )
1574 if ( nNext[i]->GetPosition()->GetDim() != 2 &&
1575 nNext[i]->GetID() < nodeOnEdge->GetID() )
1577 // look for an in-FACE node
1578 for ( int iN = 0; iN < nbN; ++iN )
1580 if ( iN == iOnE || iN == iNext[i] )
1582 SMESH_TNodeXYZ pInFace = face->GetNode( iN );
1583 gp_XYZ v = pInFace - segEnd;
1586 segVec = SMESH_TNodeXYZ( nNext[i] ) - segEnd;
1587 segLen = segVec.Modulus();
1589 double distToSeg = v.Crossed( segVec ).Modulus() / segLen;
1590 faceSize = Min( faceSize, distToSeg );
1598 //================================================================================
1600 * \brief Return direction of axis or revolution of a surface
1602 //================================================================================
1604 bool getRovolutionAxis( const Adaptor3d_Surface& surface,
1607 switch ( surface.GetType() ) {
1610 gp_Cone cone = surface.Cone();
1611 axis = cone.Axis().Direction();
1614 case GeomAbs_Sphere:
1616 gp_Sphere sphere = surface.Sphere();
1617 axis = sphere.Position().Direction();
1620 case GeomAbs_SurfaceOfRevolution:
1622 axis = surface.AxeOfRevolution().Direction();
1625 //case GeomAbs_SurfaceOfExtrusion:
1626 case GeomAbs_OffsetSurface:
1628 Handle(Adaptor3d_HSurface) base = surface.BasisSurface();
1629 return getRovolutionAxis( base->Surface(), axis );
1631 default: return false;
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
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;
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;
1665 ~PyDump() { Finish(); cout << "NB FUNCTIONS: " << theNbPyFunc << endl; }
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; }
1688 struct PyDump { PyDump(SMESH_Mesh&) {} void Finish() {} };
1689 #define dumpFunction(f) f
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 ) {}
1700 using namespace VISCOUS_3D;
1702 //================================================================================
1704 * \brief Constructor of _ViscousBuilder
1706 //================================================================================
1708 _ViscousBuilder::_ViscousBuilder()
1710 _error = SMESH_ComputeError::New(COMPERR_OK);
1714 //================================================================================
1716 * \brief Stores error description and returns false
1718 //================================================================================
1720 bool _ViscousBuilder::error(const string& text, int solidId )
1722 const string prefix = string("Viscous layers builder: ");
1723 _error->myName = COMPERR_ALGO_FAILED;
1724 _error->myComment = prefix + text;
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 )
1732 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1733 if ( smError && smError->myAlgo )
1734 _error->myAlgo = smError->myAlgo;
1736 sm->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
1738 // set KO to all solids
1739 for ( size_t i = 0; i < _sdVec.size(); ++i )
1741 if ( _sdVec[i]._index == solidId )
1743 sm = _mesh->GetSubMesh( _sdVec[i]._solid );
1744 if ( !sm->IsEmpty() )
1746 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1747 if ( !smError || smError->IsOK() )
1749 smError = SMESH_ComputeError::New( COMPERR_ALGO_FAILED, prefix + "failed");
1750 sm->ComputeStateEngine( SMESH_subMesh::CHECK_COMPUTE_STATE );
1754 makeGroupOfLE(); // debug
1759 //================================================================================
1761 * \brief At study restoration, restore event listeners used to clear an inferior
1762 * dim sub-mesh modified by viscous layers
1764 //================================================================================
1766 void _ViscousBuilder::RestoreListeners()
1771 //================================================================================
1773 * \brief computes SMESH_ProxyMesh::SubMesh::_n2n
1775 //================================================================================
1777 bool _ViscousBuilder::MakeN2NMap( _MeshOfSolid* pm )
1779 SMESH_subMesh* solidSM = pm->mySubMeshes.front();
1780 TopExp_Explorer fExp( solidSM->GetSubShape(), TopAbs_FACE );
1781 for ( ; fExp.More(); fExp.Next() )
1783 SMESHDS_SubMesh* srcSmDS = pm->GetMeshDS()->MeshElements( fExp.Current() );
1784 const SMESH_ProxyMesh::SubMesh* prxSmDS = pm->GetProxySubMesh( fExp.Current() );
1786 if ( !srcSmDS || !prxSmDS || !srcSmDS->NbElements() || !prxSmDS->NbElements() )
1788 if ( srcSmDS->GetElements()->next() == prxSmDS->GetElements()->next())
1791 if ( srcSmDS->NbElements() != prxSmDS->NbElements() )
1792 return error( "Different nb elements in a source and a proxy sub-mesh", solidSM->GetId());
1794 SMDS_ElemIteratorPtr srcIt = srcSmDS->GetElements();
1795 SMDS_ElemIteratorPtr prxIt = prxSmDS->GetElements();
1796 while( prxIt->more() )
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 );
1806 pm->_n2nMapComputed = true;
1810 //================================================================================
1812 * \brief Does its job
1814 //================================================================================
1816 SMESH_ComputeErrorPtr _ViscousBuilder::Compute(SMESH_Mesh& theMesh,
1817 const TopoDS_Shape& theShape)
1819 // TODO: set priority of solids during Gen::Compute()
1823 // check if proxy mesh already computed
1824 TopExp_Explorer exp( theShape, TopAbs_SOLID );
1826 return error("No SOLID's in theShape"), _error;
1828 if ( _ViscousListener::GetSolidMesh( _mesh, exp.Current(), /*toCreate=*/false))
1829 return SMESH_ComputeErrorPtr(); // everything already computed
1831 PyDump debugDump( theMesh );
1833 // TODO: ignore already computed SOLIDs
1834 if ( !findSolidsWithLayers())
1837 if ( !findFacesWithLayers() )
1840 for ( size_t i = 0; i < _sdVec.size(); ++i )
1842 if ( ! makeLayer(_sdVec[i]) )
1845 if ( _sdVec[i]._n2eMap.size() == 0 )
1848 if ( ! inflate(_sdVec[i]) )
1851 if ( ! refine(_sdVec[i]) )
1857 addBoundaryElements();
1859 makeGroupOfLE(); // debug
1865 //================================================================================
1867 * \brief Check validity of hypotheses
1869 //================================================================================
1871 SMESH_ComputeErrorPtr _ViscousBuilder::CheckHypotheses( SMESH_Mesh& mesh,
1872 const TopoDS_Shape& shape )
1876 if ( _ViscousListener::GetSolidMesh( _mesh, shape, /*toCreate=*/false))
1877 return SMESH_ComputeErrorPtr(); // everything already computed
1880 findSolidsWithLayers();
1881 bool ok = findFacesWithLayers( true );
1883 // remove _MeshOfSolid's of _SolidData's
1884 for ( size_t i = 0; i < _sdVec.size(); ++i )
1885 _ViscousListener::RemoveSolidMesh( _mesh, _sdVec[i]._solid );
1890 return SMESH_ComputeErrorPtr();
1893 //================================================================================
1895 * \brief Finds SOLIDs to compute using viscous layers. Fills _sdVec
1897 //================================================================================
1899 bool _ViscousBuilder::findSolidsWithLayers()
1902 TopTools_IndexedMapOfShape allSolids;
1903 TopExp::MapShapes( _mesh->GetShapeToMesh(), TopAbs_SOLID, allSolids );
1904 _sdVec.reserve( allSolids.Extent());
1906 SMESH_Gen* gen = _mesh->GetGen();
1907 SMESH_HypoFilter filter;
1908 for ( int i = 1; i <= allSolids.Extent(); ++i )
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 )))
1922 TopoDS_Shape hypShape;
1923 filter.Init( filter.Is( viscHyp ));
1924 _mesh->GetHypothesis( allSolids(i), filter, true, &hypShape );
1928 _MeshOfSolid* proxyMesh = _ViscousListener::GetSolidMesh( _mesh,
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) );
1937 soData->_hyps.push_back( viscHyp );
1938 soData->_hypShapes.push_back( hypShape );
1941 if ( _sdVec.empty() )
1943 ( SMESH_Comment(StdMeshers_ViscousLayers::GetHypType()) << " hypothesis not found",0);
1948 //================================================================================
1952 //================================================================================
1954 bool _ViscousBuilder::findFacesWithLayers(const bool onlyWith)
1956 SMESH_MesherHelper helper( *_mesh );
1957 TopExp_Explorer exp;
1958 TopTools_IndexedMapOfShape solids;
1960 // collect all faces-to-ignore defined by hyp
1961 for ( size_t i = 0; i < _sdVec.size(); ++i )
1963 solids.Add( _sdVec[i]._solid );
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 )
1973 ignoreFacesOfHyps.push_back( TFacesOfHyp( set<TGeomID>(), *hyp ));
1974 getIgnoreFaces( _sdVec[i]._solid, *hyp, *hypShape, ignoreFacesOfHyps.back().first );
1977 // fill _SolidData::_face2hyp and check compatibility of hypotheses
1978 const int nbHyps = _sdVec[i]._hyps.size();
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() )
1985 const TGeomID faceID = getMeshDS()->ShapeToIndex( exp.Current() );
1987 igFacesOfHyp = ignoreFacesOfHyps.begin();
1988 for ( ; igFacesOfHyp != ignoreFacesOfHyps.end(); ++igFacesOfHyp )
1989 if ( ! igFacesOfHyp->first.count( faceID ))
1992 return error(SMESH_Comment("Several hypotheses define "
1993 "Viscous Layers on the face #") << faceID );
1994 hyp = igFacesOfHyp->second;
1997 _sdVec[i]._face2hyp.insert( make_pair( faceID, hyp ));
1999 _sdVec[i]._ignoreFaceIds.insert( faceID );
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 )
2008 nbLayersSet.insert( igFacesOfHyp->second->GetNumberLayers() );
2010 if ( nbLayersSet.size() > 1 )
2012 for ( exp.Init( _sdVec[i]._solid, TopAbs_EDGE ); exp.More(); exp.Next() )
2014 PShapeIteratorPtr fIt = helper.GetAncestors( exp.Current(), *_mesh, TopAbs_FACE );
2015 THyp hyp1 = 0, hyp2 = 0;
2016 while( const TopoDS_Shape* face = fIt->next() )
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() )
2022 ( hyp1 ? hyp2 : hyp1 ) = f2h->second;
2025 if ( hyp1 && hyp2 &&
2026 hyp1->GetNumberLayers() != hyp2->GetNumberLayers() )
2028 return error("Two hypotheses define different number of "
2029 "viscous layers on adjacent faces");
2033 } // if ( nbHyps > 1 )
2036 _sdVec[i]._ignoreFaceIds.swap( ignoreFacesOfHyps.back().first );
2040 if ( onlyWith ) // is called to check hypotheses compatibility only
2043 // fill _SolidData::_reversedFaceIds
2044 for ( size_t i = 0; i < _sdVec.size(); ++i )
2046 exp.Init( _sdVec[i]._solid.Oriented( TopAbs_FORWARD ), TopAbs_FACE );
2047 for ( ; exp.More(); exp.Next() )
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 ))
2055 _sdVec[i]._reversedFaceIds.insert( faceID );
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 )
2065 TopExp::MapShapes(_sdVec[i]._solid, TopAbs_EDGE, shapes);
2066 for ( int iE = 1; iE <= shapes.Extent(); ++iE )
2068 const TopoDS_Shape& edge = shapes(iE);
2069 // find 2 faces sharing an edge
2071 PShapeIteratorPtr fIt = helper.GetAncestors(edge, *_mesh, TopAbs_FACE);
2072 while ( fIt->more())
2074 const TopoDS_Shape* f = fIt->next();
2075 if ( helper.IsSubShape( *f, _sdVec[i]._solid))
2076 FF[ int( !FF[0].IsNull()) ] = *f;
2078 if( FF[1].IsNull() ) continue; // seam edge can be shared by 1 FACE only
2079 // check presence of layers on them
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 ))
2092 int iSolid = solids.FindIndex( *solid );
2093 int iFace = getMeshDS()->ShapeToIndex( fWOL );
2094 if ( iSolid > 0 && !_sdVec[ iSolid-1 ]._ignoreFaceIds.count( iFace ))
2096 //_sdVec[i]._noShrinkShapes.insert( iFace );
2102 if ( !fWOL.IsNull())
2104 TGeomID edgeInd = getMeshDS()->ShapeToIndex( edge );
2105 _sdVec[i]._shrinkShape2Shape.insert( make_pair( edgeInd, fWOL ));
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 );
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 )
2120 map< TGeomID, TopoDS_Shape >::iterator e2f = _sdVec[i]._shrinkShape2Shape.begin();
2121 for ( ; e2f != _sdVec[i]._shrinkShape2Shape.end(); ++e2f )
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() )
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;
2135 for ( ; iSolid < _sdVec.size(); ++iSolid )
2137 if ( _sdVec[iSolid]._solid.IsSame( *solid ) ) {
2138 if ( _sdVec[iSolid]._shrinkShape2Shape.count( edgeID ))
2139 notShrinkFace = false;
2143 if ( notShrinkFace )
2145 _sdVec[i]._noShrinkShapes.insert( edgeID );
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() ));
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
2156 TopExp::MapShapes( fWOL, TopAbs_EDGE, shapes);
2157 for ( int iE = 1; iE <= shapes.Extent(); ++iE )
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 ))
2165 for ( int is1st = 0; is1st < 2; ++is1st )
2167 TopoDS_Vertex V = helper.IthVertex( is1st, E );
2168 if ( _sdVec[i]._noShrinkShapes.count( getMeshDS()->ShapeToIndex( V ) ))
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);
2181 } // while ( soIt->more() )
2182 } // loop on _sdVec[i]._shrinkShape2Shape
2183 } // loop on _sdVec to fill in _SolidData::_noShrinkShapes
2185 // Find the SHAPE along which to inflate _LayerEdge based on VERTEX
2187 for ( size_t i = 0; i < _sdVec.size(); ++i )
2190 TopExp::MapShapes(_sdVec[i]._solid, TopAbs_VERTEX, shapes);
2191 for ( int iV = 1; iV <= shapes.Extent(); ++iV )
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())
2200 const TopoDS_Shape* f = fIt->next();
2201 if ( helper.IsSubShape( *f, _sdVec[i]._solid ) )
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 );
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() )
2217 helper.SetSubShape( facesWOL[0] );
2218 if ( helper.IsRealSeam( vInd )) // inflate along a seam edge?
2220 TopoDS_Shape seamEdge;
2221 PShapeIteratorPtr eIt = helper.GetAncestors(vertex, *_mesh, TopAbs_EDGE);
2222 while ( eIt->more() && seamEdge.IsNull() )
2224 const TopoDS_Shape* e = eIt->next();
2225 if ( helper.IsRealSeam( *e ) )
2228 if ( !seamEdge.IsNull() )
2230 _sdVec[i]._shrinkShape2Shape.insert( make_pair( vInd, seamEdge ));
2234 _sdVec[i]._shrinkShape2Shape.insert( make_pair( vInd, facesWOL[0] ));
2239 // find an edge shared by 2 faces
2240 PShapeIteratorPtr eIt = helper.GetAncestors(vertex, *_mesh, TopAbs_EDGE);
2241 while ( eIt->more())
2243 const TopoDS_Shape* e = eIt->next();
2244 if ( helper.IsSubShape( *e, facesWOL[0]) &&
2245 helper.IsSubShape( *e, facesWOL[1]))
2247 _sdVec[i]._shrinkShape2Shape.insert( make_pair( vInd, *e )); break;
2253 return error("Not yet supported case", _sdVec[i]._index);
2258 // add FACEs of other SOLIDs to _ignoreFaceIds
2259 for ( size_t i = 0; i < _sdVec.size(); ++i )
2262 TopExp::MapShapes(_sdVec[i]._solid, TopAbs_FACE, shapes);
2264 for ( exp.Init( _mesh->GetShapeToMesh(), TopAbs_FACE ); exp.More(); exp.Next() )
2266 if ( !shapes.Contains( exp.Current() ))
2267 _sdVec[i]._ignoreFaceIds.insert( getMeshDS()->ShapeToIndex( exp.Current() ));
2274 //================================================================================
2276 * \brief Finds FACEs w/o layers for a given SOLID by an hypothesis
2278 //================================================================================
2280 void _ViscousBuilder::getIgnoreFaces(const TopoDS_Shape& solid,
2281 const StdMeshers_ViscousLayers* hyp,
2282 const TopoDS_Shape& hypShape,
2283 set<TGeomID>& ignoreFaceIds)
2285 TopExp_Explorer exp;
2287 vector<TGeomID> ids = hyp->GetBndShapes();
2288 if ( hyp->IsToIgnoreShapes() ) // FACEs to ignore are given
2290 for ( size_t ii = 0; ii < ids.size(); ++ii )
2292 const TopoDS_Shape& s = getMeshDS()->IndexToShape( ids[ii] );
2293 if ( !s.IsNull() && s.ShapeType() == TopAbs_FACE )
2294 ignoreFaceIds.insert( ids[ii] );
2297 else // FACEs with layers are given
2299 exp.Init( solid, TopAbs_FACE );
2300 for ( ; exp.More(); exp.Next() )
2302 TGeomID faceInd = getMeshDS()->ShapeToIndex( exp.Current() );
2303 if ( find( ids.begin(), ids.end(), faceInd ) == ids.end() )
2304 ignoreFaceIds.insert( faceInd );
2308 // ignore internal FACEs if inlets and outlets are specified
2309 if ( hyp->IsToIgnoreShapes() )
2311 TopTools_IndexedDataMapOfShapeListOfShape solidsOfFace;
2312 TopExp::MapShapesAndAncestors( hypShape,
2313 TopAbs_FACE, TopAbs_SOLID, solidsOfFace);
2315 for ( exp.Init( solid, TopAbs_FACE ); exp.More(); exp.Next() )
2317 const TopoDS_Face& face = TopoDS::Face( exp.Current() );
2318 if ( SMESH_MesherHelper::NbAncestors( face, *_mesh, TopAbs_SOLID ) < 2 )
2321 int nbSolids = solidsOfFace.FindFromKey( face ).Extent();
2323 ignoreFaceIds.insert( getMeshDS()->ShapeToIndex( face ));
2328 //================================================================================
2330 * \brief Create the inner surface of the viscous layer and prepare data for infation
2332 //================================================================================
2334 bool _ViscousBuilder::makeLayer(_SolidData& data)
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() )
2342 SMESH_subMesh* fSubM = _mesh->GetSubMesh( exp.Current() );
2343 if ( ! data._ignoreFaceIds.count( fSubM->GetId() ))
2344 faceIds.insert( fSubM->GetId() );
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 )
2352 TGeomID shapeInd = s2s->first;
2353 for ( size_t i = 0; i < _sdVec.size(); ++i )
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() )
2360 data._s2neMap.insert( make_pair( shapeInd, &_sdVec[i]._n2eMap ));
2366 // Create temporary faces and _LayerEdge's
2368 dumpFunction(SMESH_Comment("makeLayers_")<<data._index);
2370 data._stepSize = Precision::Infinite();
2371 data._stepSizeNodes[0] = 0;
2373 SMESH_MesherHelper helper( *_mesh );
2374 helper.SetSubShape( data._solid );
2375 helper.SetElementsOnShape( true );
2377 vector< const SMDS_MeshNode*> newNodes; // of a mesh face
2378 TNode2Edge::iterator n2e2;
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 );
2385 // set data of _EdgesOnShape's
2386 if ( SMESH_subMesh* sm = _mesh->GetSubMesh( data._solid ))
2388 SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(/*includeSelf=*/false);
2389 while ( smIt->more() )
2392 if ( sm->GetSubShape().ShapeType() == TopAbs_FACE &&
2393 !faceIds.count( sm->GetId() ))
2395 setShapeData( edgesByGeom[ sm->GetId() ], sm, data );
2398 // make _LayerEdge's
2399 for ( set<TGeomID>::iterator id = faceIds.begin(); id != faceIds.end(); ++id )
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);
2406 SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
2407 if ( !smDS ) return error(SMESH_Comment("Not meshed face ") << *id, data._index );
2409 SMDS_ElemIteratorPtr eIt = smDS->GetElements();
2410 while ( eIt->more() )
2412 const SMDS_MeshElement* face = eIt->next();
2413 double faceMaxCosin = -1;
2414 _LayerEdge* maxCosinEdge = 0;
2415 int nbDegenNodes = 0;
2417 newNodes.resize( face->NbCornerNodes() );
2418 for ( size_t i = 0 ; i < newNodes.size(); ++i )
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 );
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() )) {
2441 TNode2Edge::iterator n2e = data._n2eMap.insert( make_pair( n, (_LayerEdge*)0 )).first;
2442 if ( !(*n2e).second )
2445 _LayerEdge* edge = new _LayerEdge();
2446 edge->_nodes.push_back( n );
2448 edgesByGeom[ shapeID ]._edges.push_back( edge );
2449 const bool noShrink = data._noShrinkShapes.count( shapeID );
2451 SMESH_TNodeXYZ xyz( n );
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() ))
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() );
2470 edge->_nodes.push_back( helper.AddNode( xyz.X(), xyz.Y(), xyz.Z() ));
2472 if ( !setEdgeData( *edge, edgesByGeom[ shapeID ], helper, data ))
2475 if ( edge->_nodes.size() < 2 )
2476 edge->Block( data );
2477 //data._noShrinkShapes.insert( shapeID );
2479 dumpMove(edge->_nodes.back());
2481 if ( edge->_cosin > faceMaxCosin )
2483 faceMaxCosin = edge->_cosin;
2484 maxCosinEdge = edge;
2487 newNodes[ i ] = n2e->second->_nodes.back();
2490 data._n2eMap.insert( make_pair( face->GetNode( i ), n2e->second ));
2492 if ( newNodes.size() - nbDegenNodes < 2 )
2495 // create a temporary face
2496 const SMDS_MeshElement* newFace =
2497 new _TmpMeshFace( newNodes, --_tmpFaceID, face->getshapeId(), face->getIdInShape() );
2498 proxySub->AddElement( newFace );
2500 // compute inflation step size by min size of element on a convex surface
2501 if ( faceMaxCosin > theMinSmoothCosin )
2502 limitStepSize( data, face, maxCosinEdge );
2504 } // loop on 2D elements on a FACE
2505 } // loop on FACEs of a SOLID to create _LayerEdge's
2508 // Set _LayerEdge::_neibors
2509 TNode2Edge::iterator n2e;
2510 for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
2512 _EdgesOnShape& eos = data._edgesOnShape[iS];
2513 for ( size_t i = 0; i < eos._edges.size(); ++i )
2515 _LayerEdge* edge = eos._edges[i];
2516 TIDSortedNodeSet nearNodes;
2517 SMDS_ElemIteratorPtr fIt = edge->_nodes[0]->GetInverseElementIterator(SMDSAbs_Face);
2518 while ( fIt->more() )
2520 const SMDS_MeshElement* f = fIt->next();
2521 if ( !data._ignoreFaceIds.count( f->getshapeId() ))
2522 nearNodes.insert( f->begin_nodes(), f->end_nodes() );
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 );
2533 data._epsilon = 1e-7;
2534 if ( data._stepSize < 1. )
2535 data._epsilon *= data._stepSize;
2537 if ( !findShapesToSmooth( data ))
2540 // limit data._stepSize depending on surface curvature and fill data._convexFaces
2541 limitStepSizeByCurvature( data ); // !!! it must be before node substitution in _Simplex
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 )
2547 _EdgesOnShape& eos = data._edgesOnShape[iS];
2548 for ( size_t i = 0; i < eos._edges.size(); ++i )
2550 _LayerEdge* edge = eos._edges[i];
2551 if ( edge->IsOnEdge() )
2553 // get neighbor nodes
2554 bool hasData = ( edge->_2neibors->_edges[0] );
2555 if ( hasData ) // _LayerEdge is a copy of another one
2557 nn[0] = edge->_2neibors->srcNode(0);
2558 nn[1] = edge->_2neibors->srcNode(1);
2560 else if ( !findNeiborsOnEdge( edge, nn[0],nn[1], eos, data ))
2564 // set neighbor _LayerEdge's
2565 for ( int j = 0; j < 2; ++j )
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;
2572 edge->SetDataByNeighbors( nn[0], nn[1], eos, helper );
2575 for ( size_t j = 0; j < edge->_simplices.size(); ++j )
2577 _Simplex& s = edge->_simplices[j];
2578 s._nNext = data._n2eMap[ s._nNext ]->_nodes.back();
2579 s._nPrev = data._n2eMap[ s._nPrev ]->_nodes.back();
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() ))
2587 // Generally we should not get here
2588 if ( eos.ShapeType() != TopAbs_EDGE )
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() )
2594 const _LayerEdge* vEdge = n2e->second;
2595 edge->_normal = vEdge->_normal;
2596 edge->_lenFactor = vEdge->_lenFactor;
2597 edge->_cosin = vEdge->_cosin;
2600 } // loop on data._edgesOnShape._edges
2601 } // loop on data._edgesOnShape
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() )
2608 // if ( _EdgesOnShape* eos = data.GetShapeEdges( e2c->first ))
2609 // data.Sort2NeiborsOnEdge( eos->_edges );
2616 //================================================================================
2618 * \brief Compute inflation step size by min size of element on a convex surface
2620 //================================================================================
2622 void _ViscousBuilder::limitStepSize( _SolidData& data,
2623 const SMDS_MeshElement* face,
2624 const _LayerEdge* maxCosinEdge )
2627 double minSize = 10 * data._stepSize;
2628 const int nbNodes = face->NbCornerNodes();
2629 for ( int i = 0; i < nbNodes; ++i )
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 )
2636 double dist = SMESH_TNodeXYZ( curN ).Distance( nextN );
2637 if ( dist < minSize )
2638 minSize = dist, iN = i;
2641 double newStep = 0.8 * minSize / maxCosinEdge->_lenFactor;
2642 if ( newStep < data._stepSize )
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 ));
2651 //================================================================================
2653 * \brief Compute inflation step size by min size of element on a convex surface
2655 //================================================================================
2657 void _ViscousBuilder::limitStepSize( _SolidData& data, const double minSize )
2659 if ( minSize < data._stepSize )
2661 data._stepSize = minSize;
2662 if ( data._stepSizeNodes[0] )
2665 SMESH_TNodeXYZ(data._stepSizeNodes[0]).Distance(data._stepSizeNodes[1]);
2666 data._stepSizeCoeff = data._stepSize / dist;
2671 //================================================================================
2673 * \brief Limit data._stepSize by evaluating curvature of shapes and fill data._convexFaces
2675 //================================================================================
2677 void _ViscousBuilder::limitStepSizeByCurvature( _SolidData& data )
2679 const int nbTestPnt = 5; // on a FACE sub-shape
2681 BRepLProp_SLProps surfProp( 2, 1e-6 );
2682 SMESH_MesherHelper helper( *_mesh );
2684 data._convexFaces.clear();
2686 for ( size_t iS = 0; iS < data._edgesOnShape.size(); ++iS )
2688 _EdgesOnShape& eof = data._edgesOnShape[iS];
2689 if ( eof.ShapeType() != TopAbs_FACE ||
2690 data._ignoreFaceIds.count( eof._shapeID ))
2693 TopoDS_Face F = TopoDS::Face( eof._shape );
2694 SMESH_subMesh * sm = eof._subMesh;
2695 const TGeomID faceID = eof._shapeID;
2697 BRepAdaptor_Surface surface( F, false );
2698 surfProp.SetSurface( surface );
2700 bool isTooCurved = false;
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() )
2708 const TGeomID subID = sm->GetId();
2709 // find _LayerEdge's of a sub-shape
2711 if (( eos = data.GetShapeEdges( subID )))
2712 cnvFace._subIdToEOS.insert( make_pair( subID, eos ));
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 )
2721 gp_XY uv = helper.GetNodeUV( F, eos->_edges[ i ]->_nodes[0] );
2722 surfProp.SetParameters( uv.X(), uv.Y() );
2723 if ( !surfProp.IsCurvatureDefined() )
2725 if ( surfProp.MaxCurvature() * oriFactor > minCurvature )
2727 limitStepSize( data, 0.9 / surfProp.MaxCurvature() * oriFactor );
2730 if ( surfProp.MinCurvature() * oriFactor > minCurvature )
2732 limitStepSize( data, 0.9 / surfProp.MinCurvature() * oriFactor );
2736 } // loop on sub-shapes of the FACE
2738 if ( !isTooCurved ) continue;
2740 _ConvexFace & convFace =
2741 data._convexFaces.insert( make_pair( faceID, cnvFace )).first->second;
2744 convFace._normalsFixed = false;
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() )
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 )
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 )
2760 convFace._simplexTestEdges.push_back( ledge );
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.
2772 set< const SMDS_MeshNode* > usedNodes;
2774 // look for _LayerEdge's with null _sWOL
2775 id2eos = convFace._subIdToEOS.begin();
2776 for ( ; id2eos != convFace._subIdToEOS.end(); ++id2eos )
2778 _EdgesOnShape& eos = * id2eos->second;
2779 if ( !eos._sWOL.IsNull() )
2781 for ( size_t i = 0; i < eos._edges.size(); ++i )
2783 _LayerEdge* ledge = eos._edges[ i ];
2784 const SMDS_MeshNode* srcNode = ledge->_nodes[0];
2785 if ( !usedNodes.insert( srcNode ).second ) continue;
2787 for ( size_t i = 0; i < ledge->_simplices.size(); ++i )
2789 usedNodes.insert( ledge->_simplices[i]._nPrev );
2790 usedNodes.insert( ledge->_simplices[i]._nNext );
2792 convFace._simplexTestEdges.push_back( ledge );
2796 } // loop on FACEs of data._solid
2799 //================================================================================
2801 * \brief Detect shapes (and _LayerEdge's on them) to smooth
2803 //================================================================================
2805 bool _ViscousBuilder::findShapesToSmooth( _SolidData& data )
2807 // define allowed thickness
2808 computeGeomSize( data ); // compute data._geomSize and _LayerEdge::_maxLen
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 )
2815 data._maxThickness = Max( data._maxThickness, (*hyp)->GetTotalThickness() );
2816 data._minThickness = Min( data._minThickness, (*hyp)->GetTotalThickness() );
2818 //const double tgtThick = /*Min( 0.5 * data._geomSize, */data._maxThickness;
2820 // Find shapes needing smoothing; such a shape has _LayerEdge._normal on it's
2821 // boundry inclined to the shape at a sharp angle
2823 //list< TGeomID > shapesToSmooth;
2824 TopTools_MapOfShape edgesOfSmooFaces;
2826 SMESH_MesherHelper helper( *_mesh );
2829 vector< _EdgesOnShape >& edgesByGeom = data._edgesOnShape;
2830 data._nbShapesToSmooth = 0;
2832 for ( size_t iS = 0; iS < edgesByGeom.size(); ++iS ) // check FACEs
2834 _EdgesOnShape& eos = edgesByGeom[iS];
2835 eos._toSmooth = false;
2836 if ( eos._edges.empty() || eos.ShapeType() != TopAbs_FACE )
2839 double tgtThick = eos._hyp.GetTotalThickness();
2840 TopExp_Explorer eExp( edgesByGeom[iS]._shape, TopAbs_EDGE );
2841 for ( ; eExp.More() && !eos._toSmooth; eExp.Next() )
2843 TGeomID iE = getMeshDS()->ShapeToIndex( eExp.Current() );
2844 vector<_LayerEdge*>& eE = edgesByGeom[ iE ]._edges;
2845 if ( eE.empty() ) continue;
2848 for ( size_t i = 0; i < eE.size() && !eos._toSmooth; ++i )
2849 if ( eE[i]->_cosin > theMinSmoothCosin )
2851 SMDS_ElemIteratorPtr fIt = eE[i]->_nodes[0]->GetInverseElementIterator(SMDSAbs_Face);
2852 while ( fIt->more() && !eos._toSmooth )
2854 const SMDS_MeshElement* face = fIt->next();
2855 if ( face->getshapeId() == eos._shapeID &&
2856 getDistFromEdge( face, eE[i]->_nodes[0], faceSize ))
2858 eos._toSmooth = needSmoothing( eE[i]->_cosin, tgtThick, faceSize );
2863 if ( eos._toSmooth )
2865 for ( eExp.ReInit(); eExp.More(); eExp.Next() )
2866 edgesOfSmooFaces.Add( eExp.Current() );
2868 data.PrepareEdgesToSmoothOnFace( &edgesByGeom[iS], /*substituteSrcNodes=*/false );
2870 data._nbShapesToSmooth += eos._toSmooth;
2874 for ( size_t iS = 0; iS < edgesByGeom.size(); ++iS ) // check EDGEs
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;
2881 const TopoDS_Edge& E = TopoDS::Edge( edgesByGeom[iS]._shape );
2882 if ( SMESH_Algo::isDegenerated( E ) || !edgesOfSmooFaces.Contains( E ))
2885 double tgtThick = eos._hyp.GetTotalThickness();
2886 for ( TopoDS_Iterator vIt( E ); vIt.More() && !eos._toSmooth; vIt.Next() )
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 )
2897 // always smooth analytic EDGEs
2898 Handle(Geom_Curve) curve = _Smoother1D::CurveForSmooth( E, eos, helper );
2899 eos._toSmooth = ! curve.IsNull();
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 )
2905 const SMDS_MeshElement* endSeg = eIt->next();
2906 if ( endSeg->getshapeId() == (int) iS )
2909 SMESH_TNodeXYZ( endSeg->GetNode(0) ).Distance( endSeg->GetNode(1 ));
2910 eos._toSmooth = needSmoothing( cosinAbs, tgtThick, segLen );
2913 if ( eos._toSmooth )
2915 eos._edgeSmoother = new _Smoother1D( curve, eos );
2917 for ( size_t i = 0; i < eos._edges.size(); ++i )
2918 eos._edges[i]->Set( _LayerEdge::TO_SMOOTH );
2922 data._nbShapesToSmooth += eos._toSmooth;
2926 // Reset _cosin if no smooth is allowed by the user
2927 for ( size_t iS = 0; iS < edgesByGeom.size(); ++iS )
2929 _EdgesOnShape& eos = edgesByGeom[iS];
2930 if ( eos._edges.empty() ) continue;
2932 if ( !eos._hyp.ToSmooth() )
2933 for ( size_t i = 0; i < eos._edges.size(); ++i )
2934 eos._edges[i]->SetCosin( 0 );
2938 // Fill _eosC1 to make that C1 FACEs and EGDEs between them to be smoothed as a whole
2940 TopTools_MapOfShape c1VV;
2942 for ( size_t iS = 0; iS < edgesByGeom.size(); ++iS ) // check FACEs
2944 _EdgesOnShape& eos = edgesByGeom[iS];
2945 if ( eos._edges.empty() ||
2946 eos.ShapeType() != TopAbs_FACE ||
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() )
2955 SMESH_subMesh* sm = smQueue.front();
2956 smQueue.pop_front();
2957 SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(/*includeSelf=*/false);
2958 while ( smIt->more() )
2961 if ( sm->GetSubShape().ShapeType() == TopAbs_VERTEX )
2962 allVV.Add( sm->GetSubShape() );
2963 if ( sm->GetSubShape().ShapeType() != TopAbs_EDGE ||
2964 !checkedEE.Add( sm->GetSubShape() ))
2967 _EdgesOnShape* eoe = data.GetShapeEdges( sm->GetId() );
2968 vector<_LayerEdge*>& eE = eoe->_edges;
2969 if ( eE.empty() || !eoe->_sWOL.IsNull() )
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 );
2978 // check that mesh faces are C1 as well
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 ))
2985 while ( fIt->more() && isC1 )
2986 isC1 = ( SMESH_MeshAlgos::FaceNormal( fIt->next(), norm2, /*normalized=*/true ) &&
2987 Abs( norm1 * norm2 ) >= ( 1. - theMinSmoothCosin ));
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() )
2996 _EdgesOnShape* eof = data.GetShapeEdges( *face );
2997 if ( !eof ) continue; // other solid
2998 if ( !eos.HasC1( eoe ))
3000 eos._eosC1.push_back( eoe );
3001 eoe->_toSmooth = false;
3002 data.PrepareEdgesToSmoothOnFace( eoe, /*substituteSrcNodes=*/false );
3004 if ( eos._shapeID != eof->_shapeID && !eos.HasC1( eof ))
3006 eos._eosC1.push_back( eof );
3007 eof->_toSmooth = false;
3008 data.PrepareEdgesToSmoothOnFace( eof, /*substituteSrcNodes=*/false );
3009 smQueue.push_back( eof->_subMesh );
3014 if ( eos._eosC1.empty() )
3017 // check VERTEXes of C1 FACEs
3018 TopTools_MapIteratorOfMapOfShape vIt( allVV );
3019 for ( ; vIt.More(); vIt.Next() )
3021 _EdgesOnShape* eov = data.GetShapeEdges( vIt.Key() );
3022 if ( !eov || eov->_edges.empty() || !eov->_sWOL.IsNull() )
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() )
3029 _EdgesOnShape* eof = data.GetShapeEdges( *face );
3030 if ( !eof ) continue; // other solid
3031 isC1 = ( face->IsSame( eos._shape ) || eos.HasC1( eof ));
3037 eos._eosC1.push_back( eov );
3038 data.PrepareEdgesToSmoothOnFace( eov, /*substituteSrcNodes=*/false );
3039 c1VV.Add( eov->_shape );
3043 } // fill _eosC1 of FACEs
3048 vector< pair< _EdgesOnShape*, gp_XYZ > > dirOfEdges;
3050 for ( size_t iS = 0; iS < edgesByGeom.size(); ++iS ) // check VERTEXes
3052 _EdgesOnShape& eov = edgesByGeom[iS];
3053 if ( eov._edges.empty() ||
3054 eov.ShapeType() != TopAbs_VERTEX ||
3055 c1VV.Contains( eov._shape ))
3057 const TopoDS_Vertex& V = TopoDS::Vertex( eov._shape );
3059 // get directions of surrounding EDGEs
3061 PShapeIteratorPtr fIt = helper.GetAncestors( eov._shape, *_mesh, TopAbs_EDGE );
3062 while ( const TopoDS_Shape* e = fIt->next() )
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() ));
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 )
3076 double dot = dirOfEdges[i].second * dirOfEdges[j].second;
3077 bool isC1 = ( dot < - ( 1. - theMinSmoothCosin ));
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;
3089 } // fill _eosC1 of VERTEXes
3096 //================================================================================
3098 * \brief initialize data of _EdgesOnShape
3100 //================================================================================
3102 void _ViscousBuilder::setShapeData( _EdgesOnShape& eos,
3106 if ( !eos._shape.IsNull() ||
3107 sm->GetSubShape().ShapeType() == TopAbs_WIRE )
3110 SMESH_MesherHelper helper( *_mesh );
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;
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;
3126 eos._isRegularSWOL = true;
3127 if ( eos.SWOLType() == TopAbs_FACE )
3129 const TopoDS_Face& F = TopoDS::Face( eos._sWOL );
3130 Handle(ShapeAnalysis_Surface) surface = helper.GetSurface( F );
3131 eos._isRegularSWOL = ( ! surface->HasSingularities( 1e-7 ));
3135 if ( data._hyps.size() == 1 )
3137 eos._hyp = data._hyps.back();
3141 // compute average StdMeshers_ViscousLayers parameters
3142 map< TGeomID, const StdMeshers_ViscousLayers* >::iterator f2hyp;
3143 if ( eos.ShapeType() == TopAbs_FACE )
3145 if (( f2hyp = data._face2hyp.find( eos._shapeID )) != data._face2hyp.end() )
3146 eos._hyp = f2hyp->second;
3150 PShapeIteratorPtr fIt = helper.GetAncestors( eos._shape, *_mesh, TopAbs_FACE );
3151 while ( const TopoDS_Shape* face = fIt->next() )
3153 TGeomID faceID = getMeshDS()->ShapeToIndex( *face );
3154 if (( f2hyp = data._face2hyp.find( faceID )) != data._face2hyp.end() )
3155 eos._hyp.Add( f2hyp->second );
3161 if ( ! eos._hyp.UseSurfaceNormal() )
3163 if ( eos.ShapeType() == TopAbs_FACE ) // get normals to elements on a FACE
3165 SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
3166 eos._faceNormals.resize( smDS->NbElements() );
3168 SMDS_ElemIteratorPtr eIt = smDS->GetElements();
3169 for ( int iF = 0; eIt->more(); ++iF )
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 );
3176 if ( !helper.IsReversedSubMesh( TopoDS::Face( eos._shape )))
3177 for ( size_t iF = 0; iF < eos._faceNormals.size(); ++iF )
3178 eos._faceNormals[iF].Reverse();
3180 else // find EOS of adjacent FACEs
3182 PShapeIteratorPtr fIt = helper.GetAncestors( eos._shape, *_mesh, TopAbs_FACE );
3183 while ( const TopoDS_Shape* face = fIt->next() )
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;
3195 //================================================================================
3197 * \brief Returns normal of a face
3199 //================================================================================
3201 bool _EdgesOnShape::GetNormal( const SMDS_MeshElement* face, gp_Vec& norm )
3204 const _EdgesOnShape* eos = 0;
3206 if ( face->getshapeId() == _shapeID )
3212 for ( size_t iF = 0; iF < _faceEOS.size() && !eos; ++iF )
3213 if ( face->getshapeId() == _faceEOS[ iF ]->_shapeID )
3214 eos = _faceEOS[ iF ];
3218 ( ok = ( face->getIdInShape() < (int) eos->_faceNormals.size() )))
3220 norm = eos->_faceNormals[ face->getIdInShape() ];
3224 debugMsg( "_EdgesOnShape::Normal() failed for face "<<face->GetID()
3225 << " on _shape #" << _shapeID );
3231 //================================================================================
3233 * \brief Set data of _LayerEdge needed for smoothing
3235 //================================================================================
3237 bool _ViscousBuilder::setEdgeData(_LayerEdge& edge,
3239 SMESH_MesherHelper& helper,
3242 const SMDS_MeshNode* node = edge._nodes[0]; // source node
3245 edge._maxLen = Precision::Infinite();
3248 edge._curvature = 0;
3251 // --------------------------
3252 // Compute _normal and _cosin
3253 // --------------------------
3256 edge._lenFactor = 1.;
3257 edge._normal.SetCoord(0,0,0);
3258 _Simplex::GetSimplices( node, edge._simplices, data._ignoreFaceIds, &data );
3260 int totalNbFaces = 0;
3262 std::pair< TopoDS_Face, gp_XYZ > face2Norm[20];
3266 const bool onShrinkShape = !eos._sWOL.IsNull();
3267 const bool useGeometry = (( eos._hyp.UseSurfaceNormal() ) ||
3268 ( eos.ShapeType() != TopAbs_FACE /*&& !onShrinkShape*/ ));
3270 // get geom FACEs the node lies on
3271 //if ( useGeometry )
3273 set<TGeomID> faceIds;
3274 if ( eos.ShapeType() == TopAbs_FACE )
3276 faceIds.insert( eos._shapeID );
3280 SMDS_ElemIteratorPtr fIt = node->GetInverseElementIterator(SMDSAbs_Face);
3281 while ( fIt->more() )
3282 faceIds.insert( fIt->next()->getshapeId() );
3284 set<TGeomID>::iterator id = faceIds.begin();
3285 for ( ; id != faceIds.end(); ++id )
3287 const TopoDS_Shape& s = getMeshDS()->IndexToShape( *id );
3288 if ( s.IsNull() || s.ShapeType() != TopAbs_FACE || data._ignoreFaceIds.count( *id ))
3290 F = TopoDS::Face( s );
3291 face2Norm[ totalNbFaces ].first = F;