1 // Copyright (C) 2007-2016 CEA/DEN, EDF R&D, OPEN CASCADE
3 // Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License, or (at your option) any later version.
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
23 // NETGENPlugin : C++ implementation
24 // File : NETGENPlugin_Mesher.cxx
25 // Author : Michael Sazonov (OCN)
28 //=============================================================================
30 #include "NETGENPlugin_Mesher.hxx"
31 #include "NETGENPlugin_Hypothesis_2D.hxx"
32 #include "NETGENPlugin_SimpleHypothesis_3D.hxx"
34 #include <SMDS_FaceOfNodes.hxx>
35 #include <SMDS_LinearEdge.hxx>
36 #include <SMDS_MeshElement.hxx>
37 #include <SMDS_MeshNode.hxx>
38 #include <SMESHDS_Mesh.hxx>
39 #include <SMESH_Block.hxx>
40 #include <SMESH_Comment.hxx>
41 #include <SMESH_ComputeError.hxx>
42 #include <SMESH_ControlPnt.hxx>
43 #include <SMESH_File.hxx>
44 #include <SMESH_Gen_i.hxx>
45 #include <SMESH_Mesh.hxx>
46 #include <SMESH_MesherHelper.hxx>
47 #include <SMESH_subMesh.hxx>
48 #include <StdMeshers_QuadToTriaAdaptor.hxx>
49 #include <StdMeshers_ViscousLayers2D.hxx>
51 #include <SALOMEDS_Tool.hxx>
53 #include <utilities.h>
55 #include <BRepBuilderAPI_Copy.hxx>
56 #include <BRep_Tool.hxx>
57 #include <Bnd_B3d.hxx>
58 #include <NCollection_Map.hxx>
59 #include <Standard_ErrorHandler.hxx>
60 #include <Standard_ProgramError.hxx>
61 #include <TColStd_MapOfInteger.hxx>
63 #include <TopExp_Explorer.hxx>
64 #include <TopTools_DataMapIteratorOfDataMapOfShapeInteger.hxx>
65 #include <TopTools_DataMapIteratorOfDataMapOfShapeShape.hxx>
66 #include <TopTools_DataMapOfShapeInteger.hxx>
67 #include <TopTools_DataMapOfShapeShape.hxx>
68 #include <TopTools_MapOfShape.hxx>
71 // Netgen include files
75 #include <occgeom.hpp>
76 #include <meshing.hpp>
77 //#include <ngexception.hpp>
80 extern int OCCGenerateMesh (OCCGeometry&, Mesh*&, MeshingParameters&, int, int);
82 extern int OCCGenerateMesh (OCCGeometry&, Mesh*&, int, int, char*);
84 //extern void OCCSetLocalMeshSize(OCCGeometry & geom, Mesh & mesh);
85 extern MeshingParameters mparam;
86 extern volatile multithreadt multithread;
87 extern bool merge_solids;
89 // values used for occgeo.facemeshstatus
90 enum EFaceMeshStatus { FACE_NOT_TREATED = 0,
102 using namespace nglib;
106 #define nodeVec_ACCESS(index) ((SMDS_MeshNode*) nodeVec.at((index)))
108 #define nodeVec_ACCESS(index) ((SMDS_MeshNode*) nodeVec[index])
111 #define NGPOINT_COORDS(p) p(0),p(1),p(2)
114 // dump elements added to ng mesh
115 //#define DUMP_SEGMENTS
116 //#define DUMP_TRIANGLES
117 //#define DUMP_TRIANGLES_SCRIPT "/tmp/trias.py" //!< debug AddIntVerticesInSolids()
120 TopTools_IndexedMapOfShape ShapesWithLocalSize;
121 std::map<int,double> VertexId2LocalSize;
122 std::map<int,double> EdgeId2LocalSize;
123 std::map<int,double> FaceId2LocalSize;
124 std::map<int,double> SolidId2LocalSize;
126 std::vector<SMESHUtils::ControlPnt> ControlPoints;
127 std::set<int> ShapesWithControlPoints; // <-- allows calling SetLocalSize() several times w/o recomputing ControlPoints
129 //=============================================================================
133 //=============================================================================
135 NETGENPlugin_Mesher::NETGENPlugin_Mesher (SMESH_Mesh* mesh,
136 const TopoDS_Shape& aShape,
142 _fineness(NETGENPlugin_Hypothesis::GetDefaultFineness()),
143 _isViscousLayers2D(false),
152 SetDefaultParameters();
153 ShapesWithLocalSize.Clear();
154 VertexId2LocalSize.clear();
155 EdgeId2LocalSize.clear();
156 FaceId2LocalSize.clear();
157 SolidId2LocalSize.clear();
158 ControlPoints.clear();
159 ShapesWithControlPoints.clear();
162 //================================================================================
166 //================================================================================
168 NETGENPlugin_Mesher::~NETGENPlugin_Mesher()
176 //================================================================================
178 * Set pointer to NETGENPlugin_Mesher* field of the holder, that will be
179 * nullified at destruction of this
181 //================================================================================
183 void NETGENPlugin_Mesher::SetSelfPointer( NETGENPlugin_Mesher ** ptr )
194 //================================================================================
196 * \brief Initialize global NETGEN parameters with default values
198 //================================================================================
200 void NETGENPlugin_Mesher::SetDefaultParameters()
202 netgen::MeshingParameters& mparams = netgen::mparam;
203 // maximal mesh edge size
204 mparams.maxh = 0;//NETGENPlugin_Hypothesis::GetDefaultMaxSize();
206 // minimal number of segments per edge
207 mparams.segmentsperedge = NETGENPlugin_Hypothesis::GetDefaultNbSegPerEdge();
208 // rate of growth of size between elements
209 mparams.grading = NETGENPlugin_Hypothesis::GetDefaultGrowthRate();
210 // safety factor for curvatures (elements per radius)
211 mparams.curvaturesafety = NETGENPlugin_Hypothesis::GetDefaultNbSegPerRadius();
212 // create elements of second order
213 mparams.secondorder = NETGENPlugin_Hypothesis::GetDefaultSecondOrder();
214 // quad-dominated surface meshing
218 mparams.quad = NETGENPlugin_Hypothesis_2D::GetDefaultQuadAllowed();
219 _fineness = NETGENPlugin_Hypothesis::GetDefaultFineness();
220 mparams.uselocalh = NETGENPlugin_Hypothesis::GetDefaultSurfaceCurvature();
221 netgen::merge_solids = NETGENPlugin_Hypothesis::GetDefaultFuseEdges();
224 //=============================================================================
228 //=============================================================================
230 void SetLocalSize(TopoDS_Shape GeomShape, double LocalSize)
232 if ( GeomShape.IsNull() ) return;
233 TopAbs_ShapeEnum GeomType = GeomShape.ShapeType();
234 if (GeomType == TopAbs_COMPOUND) {
235 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()) {
236 SetLocalSize(it.Value(), LocalSize);
241 if (! ShapesWithLocalSize.Contains(GeomShape))
242 key = ShapesWithLocalSize.Add(GeomShape);
244 key = ShapesWithLocalSize.FindIndex(GeomShape);
245 if (GeomType == TopAbs_VERTEX) {
246 VertexId2LocalSize[key] = LocalSize;
247 } else if (GeomType == TopAbs_EDGE) {
248 EdgeId2LocalSize[key] = LocalSize;
249 } else if (GeomType == TopAbs_FACE) {
250 FaceId2LocalSize[key] = LocalSize;
251 } else if (GeomType == TopAbs_SOLID) {
252 SolidId2LocalSize[key] = LocalSize;
256 //=============================================================================
258 * Pass parameters to NETGEN
260 //=============================================================================
261 void NETGENPlugin_Mesher::SetParameters(const NETGENPlugin_Hypothesis* hyp)
265 netgen::MeshingParameters& mparams = netgen::mparam;
266 // Initialize global NETGEN parameters:
267 // maximal mesh segment size
268 mparams.maxh = hyp->GetMaxSize();
269 // maximal mesh element linear size
270 mparams.minh = hyp->GetMinSize();
271 // minimal number of segments per edge
272 mparams.segmentsperedge = hyp->GetNbSegPerEdge();
273 // rate of growth of size between elements
274 mparams.grading = hyp->GetGrowthRate();
275 // safety factor for curvatures (elements per radius)
276 mparams.curvaturesafety = hyp->GetNbSegPerRadius();
277 // create elements of second order
278 mparams.secondorder = hyp->GetSecondOrder() ? 1 : 0;
279 // quad-dominated surface meshing
280 mparams.quad = hyp->GetQuadAllowed() ? 1 : 0;
281 _optimize = hyp->GetOptimize();
282 _fineness = hyp->GetFineness();
283 mparams.uselocalh = hyp->GetSurfaceCurvature();
284 netgen::merge_solids = hyp->GetFuseEdges();
287 mparams.meshsizefilename= hyp->GetMeshSizeFile().empty() ? 0 : hyp->GetMeshSizeFile().c_str();
289 SMESH_Gen_i* smeshGen_i = SMESH_Gen_i::GetSMESHGen();
290 CORBA::Object_var anObject = smeshGen_i->GetNS()->Resolve("/myStudyManager");
291 SALOMEDS::StudyManager_var aStudyMgr = SALOMEDS::StudyManager::_narrow(anObject);
292 SALOMEDS::Study_var myStudy = aStudyMgr->GetStudyByID(hyp->GetStudyId());
293 if ( !myStudy->_is_nil() )
295 const NETGENPlugin_Hypothesis::TLocalSize localSizes = hyp->GetLocalSizesAndEntries();
296 NETGENPlugin_Hypothesis::TLocalSize::const_iterator it = localSizes.begin();
297 for ( ; it != localSizes.end() ; it++)
299 std::string entry = (*it).first;
300 double val = (*it).second;
302 GEOM::GEOM_Object_var aGeomObj;
303 SALOMEDS::SObject_var aSObj = myStudy->FindObjectID( entry.c_str() );
304 if ( !aSObj->_is_nil() ) {
305 CORBA::Object_var obj = aSObj->GetObject();
306 aGeomObj = GEOM::GEOM_Object::_narrow(obj);
309 TopoDS_Shape S = smeshGen_i->GeomObjectToShape( aGeomObj.in() );
310 ::SetLocalSize(S, val);
316 //=============================================================================
318 * Pass simple parameters to NETGEN
320 //=============================================================================
322 void NETGENPlugin_Mesher::SetParameters(const NETGENPlugin_SimpleHypothesis_2D* hyp)
326 SetDefaultParameters();
329 //=============================================================================
331 * Link - a pair of integer numbers
333 //=============================================================================
337 Link(int _n1, int _n2) : n1(_n1), n2(_n2) {}
338 Link() : n1(0), n2(0) {}
339 bool Contains( int n ) const { return n == n1 || n == n2; }
340 bool IsConnected( const Link& other ) const
342 return (( Contains( other.n1 ) || Contains( other.n2 )) && ( this != &other ));
346 int HashCode(const Link& aLink, int aLimit)
348 return HashCode(aLink.n1 + aLink.n2, aLimit);
351 Standard_Boolean IsEqual(const Link& aLink1, const Link& aLink2)
353 return (( aLink1.n1 == aLink2.n1 && aLink1.n2 == aLink2.n2 ) ||
354 ( aLink1.n1 == aLink2.n2 && aLink1.n2 == aLink2.n1 ));
359 //================================================================================
361 * \brief return id of netgen point corresponding to SMDS node
363 //================================================================================
364 typedef map< const SMDS_MeshNode*, int > TNode2IdMap;
366 int ngNodeId( const SMDS_MeshNode* node,
367 netgen::Mesh& ngMesh,
368 TNode2IdMap& nodeNgIdMap)
370 int newNgId = ngMesh.GetNP() + 1;
372 TNode2IdMap::iterator node_id = nodeNgIdMap.insert( make_pair( node, newNgId )).first;
374 if ( node_id->second == newNgId)
376 #if defined(DUMP_SEGMENTS) || defined(DUMP_TRIANGLES)
377 cout << "Ng " << newNgId << " - " << node;
379 netgen::MeshPoint p( netgen::Point<3> (node->X(), node->Y(), node->Z()) );
380 ngMesh.AddPoint( p );
382 return node_id->second;
385 //================================================================================
387 * \brief Return computed EDGEs connected to the given one
389 //================================================================================
391 list< TopoDS_Edge > getConnectedEdges( const TopoDS_Edge& edge,
392 const TopoDS_Face& face,
393 const set< SMESH_subMesh* > & computedSM,
394 const SMESH_MesherHelper& helper,
395 map< SMESH_subMesh*, set< int > >& addedEdgeSM2Faces)
398 list< TopoDS_Edge > edges;
399 list< int > nbEdgesInWire;
400 /*int nbWires =*/ SMESH_Block::GetOrderedEdges( face, edges, nbEdgesInWire);
402 // find <edge> within <edges>
403 list< TopoDS_Edge >::iterator eItFwd = edges.begin();
404 for ( ; eItFwd != edges.end(); ++eItFwd )
405 if ( edge.IsSame( *eItFwd ))
407 if ( eItFwd == edges.end()) return list< TopoDS_Edge>();
409 if ( eItFwd->Orientation() >= TopAbs_INTERNAL )
411 // connected INTERNAL edges returned from GetOrderedEdges() are wrongly oriented
412 // so treat each INTERNAL edge separately
413 TopoDS_Edge e = *eItFwd;
415 edges.push_back( e );
419 // get all computed EDGEs connected to <edge>
421 list< TopoDS_Edge >::iterator eItBack = eItFwd, ePrev;
422 TopoDS_Vertex vCommon;
423 TopTools_MapOfShape eAdded; // map used not to add a seam edge twice to <edges>
426 // put edges before <edge> to <edges> back
427 while ( edges.begin() != eItFwd )
428 edges.splice( edges.end(), edges, edges.begin() );
432 while ( ++eItFwd != edges.end() )
434 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( *eItFwd );
436 bool connected = TopExp::CommonVertex( *ePrev, *eItFwd, vCommon );
437 bool computed = sm->IsMeshComputed();
438 bool added = addedEdgeSM2Faces[ sm ].count( helper.GetSubShapeID() );
439 bool doubled = !eAdded.Add( *eItFwd );
440 bool orientOK = (( ePrev ->Orientation() < TopAbs_INTERNAL ) ==
441 ( eItFwd->Orientation() < TopAbs_INTERNAL ) );
442 if ( !connected || !computed || !orientOK || added || doubled )
444 // stop advancement; move edges from tail to head
445 while ( edges.back() != *ePrev )
446 edges.splice( edges.begin(), edges, --edges.end() );
452 while ( eItBack != edges.begin() )
456 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( *eItBack );
458 bool connected = TopExp::CommonVertex( *ePrev, *eItBack, vCommon );
459 bool computed = sm->IsMeshComputed();
460 bool added = addedEdgeSM2Faces[ sm ].count( helper.GetSubShapeID() );
461 bool doubled = !eAdded.Add( *eItBack );
462 bool orientOK = (( ePrev ->Orientation() < TopAbs_INTERNAL ) ==
463 ( eItBack->Orientation() < TopAbs_INTERNAL ) );
464 if ( !connected || !computed || !orientOK || added || doubled)
467 edges.erase( edges.begin(), ePrev );
471 if ( edges.front() != edges.back() )
473 // assure that the 1st vertex is meshed
474 TopoDS_Edge eLast = edges.back();
475 while ( !SMESH_Algo::VertexNode( SMESH_MesherHelper::IthVertex( 0, edges.front()), helper.GetMeshDS())
477 edges.front() != eLast )
478 edges.splice( edges.end(), edges, edges.begin() );
483 //================================================================================
485 * \brief Make triangulation of a shape precise enough
487 //================================================================================
489 void updateTriangulation( const TopoDS_Shape& shape )
491 // static set< Poly_Triangulation* > updated;
493 // TopLoc_Location loc;
494 // TopExp_Explorer fExp( shape, TopAbs_FACE );
495 // for ( ; fExp.More(); fExp.Next() )
497 // Handle(Poly_Triangulation) triangulation =
498 // BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
499 // if ( triangulation.IsNull() ||
500 // updated.insert( triangulation.operator->() ).second )
502 // BRepTools::Clean (shape);
505 BRepMesh_IncrementalMesh e(shape, 0.01, true);
507 catch (Standard_Failure)
510 // updated.erase( triangulation.operator->() );
511 // triangulation = BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
512 // updated.insert( triangulation.operator->() );
516 //================================================================================
518 * \brief Returns a medium node either existing in SMESH of created by NETGEN
519 * \param [in] corner1 - corner node 1
520 * \param [in] corner2 - corner node 2
521 * \param [in] defaultMedium - the node created by NETGEN
522 * \param [in] helper - holder of medium nodes existing in SMESH
523 * \return const SMDS_MeshNode* - the result node
525 //================================================================================
527 const SMDS_MeshNode* mediumNode( const SMDS_MeshNode* corner1,
528 const SMDS_MeshNode* corner2,
529 const SMDS_MeshNode* defaultMedium,
530 const SMESH_MesherHelper* helper)
534 TLinkNodeMap::const_iterator l2n =
535 helper->GetTLinkNodeMap().find( SMESH_TLink( corner1, corner2 ));
536 if ( l2n != helper->GetTLinkNodeMap().end() )
537 defaultMedium = l2n->second;
539 return defaultMedium;
542 //================================================================================
544 * \brief Assure that mesh on given shapes is quadratic
546 //================================================================================
548 // void makeQuadratic( const TopTools_IndexedMapOfShape& shapes,
549 // SMESH_Mesh* mesh )
551 // for ( int i = 1; i <= shapes.Extent(); ++i )
553 // SMESHDS_SubMesh* smDS = mesh->GetMeshDS()->MeshElements( shapes(i) );
554 // if ( !smDS ) continue;
555 // SMDS_ElemIteratorPtr elemIt = smDS->GetElements();
556 // if ( !elemIt->more() ) continue;
557 // const SMDS_MeshElement* e = elemIt->next();
558 // if ( !e || e->IsQuadratic() )
561 // TIDSortedElemSet elems;
562 // elems.insert( e );
563 // while ( elemIt->more() )
564 // elems.insert( elems.end(), elemIt->next() );
566 // SMESH_MeshEditor( mesh ).ConvertToQuadratic( /*3d=*/false, elems, /*biQuad=*/false );
570 //================================================================================
572 * \brief Restrict size of elements on the given edge
574 //================================================================================
576 void setLocalSize(const TopoDS_Edge& edge,
580 if ( size <= std::numeric_limits<double>::min() )
582 Standard_Real u1, u2;
583 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, u1, u2);
584 if ( curve.IsNull() )
586 TopoDS_Iterator vIt( edge );
587 if ( !vIt.More() ) return;
588 gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex( vIt.Value() ));
589 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
593 const int nb = (int)( 1.5 * SMESH_Algo::EdgeLength( edge ) / size );
594 Standard_Real delta = (u2-u1)/nb;
595 for(int i=0; i<nb; i++)
597 Standard_Real u = u1 + delta*i;
598 gp_Pnt p = curve->Value(u);
599 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
600 netgen::Point3d pi(p.X(), p.Y(), p.Z());
601 double resultSize = mesh.GetH(pi);
602 if ( resultSize - size > 0.1*size )
603 // netgen does restriction iff oldH/newH > 1.2 (localh.cpp:136)
604 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), resultSize/1.201 );
610 //================================================================================
612 * \brief Set local size on shapes defined by SetParameters()
614 //================================================================================
616 void NETGENPlugin_Mesher::SetLocalSize( netgen::OCCGeometry& occgeo,
617 netgen::Mesh& ngMesh )
619 for(std::map<int,double>::const_iterator it=EdgeId2LocalSize.begin(); it!=EdgeId2LocalSize.end(); it++)
621 int key = (*it).first;
622 double hi = (*it).second;
623 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
624 setLocalSize( TopoDS::Edge(shape), hi, ngMesh );
626 for(std::map<int,double>::const_iterator it=VertexId2LocalSize.begin(); it!=VertexId2LocalSize.end(); it++)
628 int key = (*it).first;
629 double hi = (*it).second;
630 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
631 gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex(shape) );
632 NETGENPlugin_Mesher::RestrictLocalSize( ngMesh, p.XYZ(), hi );
634 for(map<int,double>::const_iterator it=FaceId2LocalSize.begin(); it!=FaceId2LocalSize.end(); it++)
636 int key = (*it).first;
637 double val = (*it).second;
638 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
639 int faceNgID = occgeo.fmap.FindIndex(shape);
642 occgeo.SetFaceMaxH(faceNgID, val);
643 for ( TopExp_Explorer edgeExp( shape, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
644 setLocalSize( TopoDS::Edge( edgeExp.Current() ), val, ngMesh );
646 else if ( !ShapesWithControlPoints.count( key ))
648 SMESHUtils::createPointsSampleFromFace( TopoDS::Face( shape ), val, ControlPoints );
649 ShapesWithControlPoints.insert( key );
652 for(map<int,double>::const_iterator it=SolidId2LocalSize.begin(); it!=SolidId2LocalSize.end(); it++)
654 int key = (*it).first;
655 double val = (*it).second;
656 if ( !ShapesWithControlPoints.count( key ))
658 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
659 SMESHUtils::createPointsSampleFromSolid( TopoDS::Solid( shape ), val, ControlPoints );
660 ShapesWithControlPoints.insert( key );
664 if ( !ControlPoints.empty() )
666 for ( size_t i = 0; i < ControlPoints.size(); ++i )
667 NETGENPlugin_Mesher::RestrictLocalSize( ngMesh, ControlPoints[i].XYZ(), ControlPoints[i].Size() );
671 //================================================================================
673 * \brief Initialize netgen::OCCGeometry with OCCT shape
675 //================================================================================
677 void NETGENPlugin_Mesher::PrepareOCCgeometry(netgen::OCCGeometry& occgeo,
678 const TopoDS_Shape& shape,
680 list< SMESH_subMesh* > * meshedSM,
681 NETGENPlugin_Internals* intern)
683 updateTriangulation( shape );
686 BRepBndLib::Add (shape, bb);
687 double x1,y1,z1,x2,y2,z2;
688 bb.Get (x1,y1,z1,x2,y2,z2);
689 netgen::Point<3> p1 = netgen::Point<3> (x1,y1,z1);
690 netgen::Point<3> p2 = netgen::Point<3> (x2,y2,z2);
691 occgeo.boundingbox = netgen::Box<3> (p1,p2);
693 occgeo.shape = shape;
696 // fill maps of shapes of occgeo with not yet meshed subshapes
698 // get root submeshes
699 list< SMESH_subMesh* > rootSM;
700 const int shapeID = mesh.GetMeshDS()->ShapeToIndex( shape );
701 if ( shapeID > 0 ) { // SMESH_subMesh with ID 0 may exist, don't use it!
702 rootSM.push_back( mesh.GetSubMesh( shape ));
705 for ( TopoDS_Iterator it( shape ); it.More(); it.Next() )
706 rootSM.push_back( mesh.GetSubMesh( it.Value() ));
711 // add subshapes of empty submeshes
712 list< SMESH_subMesh* >::iterator rootIt = rootSM.begin(), rootEnd = rootSM.end();
713 for ( ; rootIt != rootEnd; ++rootIt ) {
714 SMESH_subMesh * root = *rootIt;
715 SMESH_subMeshIteratorPtr smIt = root->getDependsOnIterator(/*includeSelf=*/true,
716 /*complexShapeFirst=*/true);
717 // to find a right orientation of subshapes (PAL20462)
718 TopTools_IndexedMapOfShape subShapes;
719 TopExp::MapShapes(root->GetSubShape(), subShapes);
720 while ( smIt->more() )
722 SMESH_subMesh* sm = smIt->next();
723 TopoDS_Shape shape = sm->GetSubShape();
724 totNbFaces += ( shape.ShapeType() == TopAbs_FACE );
725 if ( intern && intern->isShapeToPrecompute( shape ))
727 if ( !meshedSM || sm->IsEmpty() )
729 if ( shape.ShapeType() != TopAbs_VERTEX )
730 shape = subShapes( subShapes.FindIndex( shape ));// shape -> index -> oriented shape
731 if ( shape.Orientation() >= TopAbs_INTERNAL )
732 shape.Orientation( TopAbs_FORWARD ); // isuue 0020676
733 switch ( shape.ShapeType() ) {
734 case TopAbs_FACE : occgeo.fmap.Add( shape ); break;
735 case TopAbs_EDGE : occgeo.emap.Add( shape ); break;
736 case TopAbs_VERTEX: occgeo.vmap.Add( shape ); break;
737 case TopAbs_SOLID :occgeo.somap.Add( shape ); break;
741 // collect submeshes of meshed shapes
744 const int dim = SMESH_Gen::GetShapeDim( shape );
745 meshedSM[ dim ].push_back( sm );
749 occgeo.facemeshstatus.SetSize (totNbFaces);
750 occgeo.facemeshstatus = 0;
751 occgeo.face_maxh_modified.SetSize(totNbFaces);
752 occgeo.face_maxh_modified = 0;
753 occgeo.face_maxh.SetSize(totNbFaces);
754 occgeo.face_maxh = netgen::mparam.maxh;
757 //================================================================================
759 * \brief Return a default min size value suitable for the given geometry.
761 //================================================================================
763 double NETGENPlugin_Mesher::GetDefaultMinSize(const TopoDS_Shape& geom,
764 const double maxSize)
766 updateTriangulation( geom );
770 const int* pi[4] = { &i1, &i2, &i3, &i1 };
773 TopExp_Explorer fExp( geom, TopAbs_FACE );
774 for ( ; fExp.More(); fExp.Next() )
776 Handle(Poly_Triangulation) triangulation =
777 BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
778 if ( triangulation.IsNull() ) continue;
779 const double fTol = BRep_Tool::Tolerance( TopoDS::Face( fExp.Current() ));
780 const TColgp_Array1OfPnt& points = triangulation->Nodes();
781 const Poly_Array1OfTriangle& trias = triangulation->Triangles();
782 for ( int iT = trias.Lower(); iT <= trias.Upper(); ++iT )
784 trias(iT).Get( i1, i2, i3 );
785 for ( int j = 0; j < 3; ++j )
787 double dist2 = points(*pi[j]).SquareDistance( points( *pi[j+1] ));
788 if ( dist2 < minh && fTol*fTol < dist2 )
790 bb.Add( points(*pi[j]));
794 if ( minh > 0.25 * bb.SquareExtent() ) // simple geometry, rough triangulation
796 minh = 1e-3 * sqrt( bb.SquareExtent());
797 //cout << "BND BOX minh = " <<minh << endl;
801 minh = 3 * sqrt( minh ); // triangulation for visualization is rather fine
802 //cout << "TRIANGULATION minh = " <<minh << endl;
804 if ( minh > 0.5 * maxSize )
810 //================================================================================
812 * \brief Restrict size of elements at a given point
814 //================================================================================
816 void NETGENPlugin_Mesher::RestrictLocalSize(netgen::Mesh& ngMesh,
819 const bool overrideMinH)
821 if ( size <= std::numeric_limits<double>::min() )
823 if ( netgen::mparam.minh > size )
827 ngMesh.SetMinimalH( size );
828 netgen::mparam.minh = size;
832 size = netgen::mparam.minh;
835 netgen::Point3d pi(p.X(), p.Y(), p.Z());
836 ngMesh.RestrictLocalH( pi, size );
839 //================================================================================
841 * \brief fill ngMesh with nodes and elements of computed submeshes
843 //================================================================================
845 bool NETGENPlugin_Mesher::FillNgMesh(netgen::OCCGeometry& occgeom,
846 netgen::Mesh& ngMesh,
847 vector<const SMDS_MeshNode*>& nodeVec,
848 const list< SMESH_subMesh* > & meshedSM,
849 SMESH_MesherHelper* quadHelper,
850 SMESH_ProxyMesh::Ptr proxyMesh)
852 TNode2IdMap nodeNgIdMap;
853 for ( size_t i = 1; i < nodeVec.size(); ++i )
854 nodeNgIdMap.insert( make_pair( nodeVec[i], i ));
856 TopTools_MapOfShape visitedShapes;
857 map< SMESH_subMesh*, set< int > > visitedEdgeSM2Faces;
858 set< SMESH_subMesh* > computedSM( meshedSM.begin(), meshedSM.end() );
860 SMESH_MesherHelper helper (*_mesh);
862 int faceNgID = ngMesh.GetNFD();
864 list< SMESH_subMesh* >::const_iterator smIt, smEnd = meshedSM.end();
865 for ( smIt = meshedSM.begin(); smIt != smEnd; ++smIt )
867 SMESH_subMesh* sm = *smIt;
868 if ( !visitedShapes.Add( sm->GetSubShape() ))
871 const SMESHDS_SubMesh * smDS = sm->GetSubMeshDS();
872 if ( !smDS ) continue;
874 switch ( sm->GetSubShape().ShapeType() )
876 case TopAbs_EDGE: { // EDGE
877 // ----------------------
878 TopoDS_Edge geomEdge = TopoDS::Edge( sm->GetSubShape() );
879 if ( geomEdge.Orientation() >= TopAbs_INTERNAL )
880 geomEdge.Orientation( TopAbs_FORWARD ); // issue 0020676
882 // Add ng segments for each not meshed FACE the EDGE bounds
883 PShapeIteratorPtr fIt = helper.GetAncestors( geomEdge, *sm->GetFather(), TopAbs_FACE );
884 while ( const TopoDS_Shape * anc = fIt->next() )
886 faceNgID = occgeom.fmap.FindIndex( *anc );
888 continue; // meshed face
890 int faceSMDSId = helper.GetMeshDS()->ShapeToIndex( *anc );
891 if ( visitedEdgeSM2Faces[ sm ].count( faceSMDSId ))
892 continue; // already treated EDGE
894 TopoDS_Face face = TopoDS::Face( occgeom.fmap( faceNgID ));
895 if ( face.Orientation() >= TopAbs_INTERNAL )
896 face.Orientation( TopAbs_FORWARD ); // issue 0020676
898 // get all meshed EDGEs of the FACE connected to geomEdge (issue 0021140)
899 helper.SetSubShape( face );
900 list< TopoDS_Edge > edges = getConnectedEdges( geomEdge, face, computedSM, helper,
901 visitedEdgeSM2Faces );
903 continue; // wrong ancestor?
905 // find out orientation of <edges> within <face>
906 TopoDS_Edge eNotSeam = edges.front();
907 if ( helper.HasSeam() )
909 list< TopoDS_Edge >::iterator eIt = edges.begin();
910 while ( helper.IsRealSeam( *eIt )) ++eIt;
911 if ( eIt != edges.end() )
914 TopAbs_Orientation fOri = helper.GetSubShapeOri( face, eNotSeam );
915 bool isForwad = ( fOri == eNotSeam.Orientation() || fOri >= TopAbs_INTERNAL );
917 // get all nodes from connected <edges>
918 const bool isQuad = smDS->IsQuadratic();
919 StdMeshers_FaceSide fSide( face, edges, _mesh, isForwad, isQuad );
920 const vector<UVPtStruct>& points = fSide.GetUVPtStruct();
921 if ( points.empty() )
922 return false; // invalid node params?
923 int i, nbSeg = fSide.NbSegments();
925 // remember EDGEs of fSide to treat only once
926 for ( int iE = 0; iE < fSide.NbEdges(); ++iE )
927 visitedEdgeSM2Faces[ helper.GetMesh()->GetSubMesh( fSide.Edge(iE )) ].insert(faceSMDSId);
929 double otherSeamParam = 0;
934 int prevNgId = ngNodeId( points[0].node, ngMesh, nodeNgIdMap );
936 for ( i = 0; i < nbSeg; ++i )
938 const UVPtStruct& p1 = points[ i ];
939 const UVPtStruct& p2 = points[ i+1 ];
941 if ( p1.node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX ) //an EDGE begins
944 if ( helper.IsRealSeam( p1.node->getshapeId() ))
946 TopoDS_Edge e = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
947 isSeam = helper.IsRealSeam( e );
950 otherSeamParam = helper.GetOtherParam( helper.GetPeriodicIndex() & 1 ? p2.u : p2.v );
957 seg[1] = prevNgId = ngNodeId( p2.node, ngMesh, nodeNgIdMap );
958 // node param on curve
959 seg.epgeominfo[ 0 ].dist = p1.param;
960 seg.epgeominfo[ 1 ].dist = p2.param;
962 seg.epgeominfo[ 0 ].u = p1.u;
963 seg.epgeominfo[ 0 ].v = p1.v;
964 seg.epgeominfo[ 1 ].u = p2.u;
965 seg.epgeominfo[ 1 ].v = p2.v;
967 //geomEdge = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
968 //seg.epgeominfo[ 0 ].edgenr = seg.epgeominfo[ 1 ].edgenr = occgeom.emap.FindIndex( geomEdge );
970 //seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
971 seg.si = faceNgID; // = geom.fmap.FindIndex (face);
972 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
973 ngMesh.AddSegment (seg);
975 SMESH_TNodeXYZ np1( p1.node ), np2( p2.node );
976 RestrictLocalSize( ngMesh, 0.5*(np1+np2), (np1-np2).Modulus() );
979 cout << "Segment: " << seg.edgenr << " on SMESH face " << helper.GetMeshDS()->ShapeToIndex( face ) << endl
980 << "\tface index: " << seg.si << endl
981 << "\tp1: " << seg[0] << endl
982 << "\tp2: " << seg[1] << endl
983 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
984 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
985 //<< "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
986 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
987 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl;
988 //<< "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
992 if ( helper.GetPeriodicIndex() && 1 ) {
993 seg.epgeominfo[ 0 ].u = otherSeamParam;
994 seg.epgeominfo[ 1 ].u = otherSeamParam;
995 swap (seg.epgeominfo[0].v, seg.epgeominfo[1].v);
997 seg.epgeominfo[ 0 ].v = otherSeamParam;
998 seg.epgeominfo[ 1 ].v = otherSeamParam;
999 swap (seg.epgeominfo[0].u, seg.epgeominfo[1].u);
1001 swap( seg[0], seg[1] );
1002 swap( seg.epgeominfo[0].dist, seg.epgeominfo[1].dist );
1003 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1004 ngMesh.AddSegment( seg );
1005 #ifdef DUMP_SEGMENTS
1006 cout << "Segment: " << seg.edgenr << endl
1007 << "\t is SEAM (reverse) of the previous. "
1008 << " Other " << (helper.GetPeriodicIndex() && 1 ? "U" : "V")
1009 << " = " << otherSeamParam << endl;
1012 else if ( fOri == TopAbs_INTERNAL )
1014 swap( seg[0], seg[1] );
1015 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1016 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1017 ngMesh.AddSegment( seg );
1018 #ifdef DUMP_SEGMENTS
1019 cout << "Segment: " << seg.edgenr << endl << "\t is REVERSE of the previous" << endl;
1023 } // loop on geomEdge ancestors
1025 if ( quadHelper ) // remember medium nodes of sub-meshes
1027 SMDS_ElemIteratorPtr edges = smDS->GetElements();
1028 while ( edges->more() )
1030 const SMDS_MeshElement* e = edges->next();
1031 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshEdge*>( e )))
1037 } // case TopAbs_EDGE
1039 case TopAbs_FACE: { // FACE
1040 // ----------------------
1041 const TopoDS_Face& geomFace = TopoDS::Face( sm->GetSubShape() );
1042 helper.SetSubShape( geomFace );
1043 bool isInternalFace = ( geomFace.Orientation() == TopAbs_INTERNAL );
1045 // Find solids the geomFace bounds
1046 int solidID1 = 0, solidID2 = 0;
1047 StdMeshers_QuadToTriaAdaptor* quadAdaptor =
1048 dynamic_cast<StdMeshers_QuadToTriaAdaptor*>( proxyMesh.get() );
1051 solidID1 = occgeom.somap.FindIndex( quadAdaptor->GetShape() );
1055 PShapeIteratorPtr solidIt = helper.GetAncestors( geomFace, *sm->GetFather(), TopAbs_SOLID);
1056 while ( const TopoDS_Shape * solid = solidIt->next() )
1058 int id = occgeom.somap.FindIndex ( *solid );
1059 if ( solidID1 && id != solidID1 ) solidID2 = id;
1063 // Add ng face descriptors of meshed faces
1065 ngMesh.AddFaceDescriptor( netgen::FaceDescriptor( faceNgID, solidID1, solidID2, 0 ));
1067 // if second oreder is required, even already meshed faces must be passed to NETGEN
1068 int fID = occgeom.fmap.Add( geomFace );
1069 if ( occgeom.facemeshstatus.Size() < fID ) occgeom.facemeshstatus.SetSize( fID );
1070 occgeom.facemeshstatus[ fID-1 ] = netgen::FACE_MESHED_OK;
1071 while ( fID < faceNgID ) // geomFace is already in occgeom.fmap, add a copy
1073 fID = occgeom.fmap.Add( BRepBuilderAPI_Copy( geomFace, /*copyGeom=*/false ));
1074 if ( occgeom.facemeshstatus.Size() < fID ) occgeom.facemeshstatus.SetSize( fID );
1075 occgeom.facemeshstatus[ fID-1 ] = netgen::FACE_MESHED_OK;
1077 // Problem with the second order in a quadrangular mesh remains.
1078 // 1) All quadrangles generated by NETGEN are moved to an inexistent face
1079 // by FillSMesh() (find "AddFaceDescriptor")
1080 // 2) Temporary triangles generated by StdMeshers_QuadToTriaAdaptor
1081 // are on faces where quadrangles were.
1082 // Due to these 2 points, wrong geom faces are used while conversion to quadratic
1083 // of the mentioned above quadrangles and triangles
1085 // Orient the face correctly in solidID1 (issue 0020206)
1086 bool reverse = false;
1088 TopoDS_Shape solid = occgeom.somap( solidID1 );
1089 TopAbs_Orientation faceOriInSolid = helper.GetSubShapeOri( solid, geomFace );
1090 if ( faceOriInSolid >= 0 )
1092 helper.IsReversedSubMesh( TopoDS::Face( geomFace.Oriented( faceOriInSolid )));
1095 // Add surface elements
1097 netgen::Element2d tri(3);
1098 tri.SetIndex( faceNgID );
1099 SMESH_TNodeXYZ xyz[3];
1101 #ifdef DUMP_TRIANGLES
1102 cout << "SMESH face " << helper.GetMeshDS()->ShapeToIndex( geomFace )
1103 << " internal="<<isInternalFace << endl;
1106 smDS = proxyMesh->GetSubMesh( geomFace );
1108 SMDS_ElemIteratorPtr faces = smDS->GetElements();
1109 while ( faces->more() )
1111 const SMDS_MeshElement* f = faces->next();
1112 if ( f->NbNodes() % 3 != 0 ) // not triangle
1114 PShapeIteratorPtr solidIt=helper.GetAncestors(geomFace,*sm->GetFather(),TopAbs_SOLID);
1115 if ( const TopoDS_Shape * solid = solidIt->next() )
1116 sm = _mesh->GetSubMesh( *solid );
1117 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1118 smError.reset( new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH,"Not triangle sub-mesh"));
1119 smError->myBadElements.push_back( f );
1123 for ( int i = 0; i < 3; ++i )
1125 const SMDS_MeshNode* node = f->GetNode( i ), * inFaceNode=0;
1128 // get node UV on face
1129 int shapeID = node->getshapeId();
1130 if ( helper.IsSeamShape( shapeID ))
1132 if ( helper.IsSeamShape( f->GetNodeWrap( i+1 )->getshapeId() ))
1133 inFaceNode = f->GetNodeWrap( i-1 );
1135 inFaceNode = f->GetNodeWrap( i+1 );
1137 gp_XY uv = helper.GetNodeUV( geomFace, node, inFaceNode );
1139 int ind = reverse ? 3-i : i+1;
1140 tri.GeomInfoPi(ind).u = uv.X();
1141 tri.GeomInfoPi(ind).v = uv.Y();
1142 tri.PNum (ind) = ngNodeId( node, ngMesh, nodeNgIdMap );
1145 // pass a triangle size to NG size-map
1146 double size = ( ( xyz[0] - xyz[1] ).Modulus() +
1147 ( xyz[1] - xyz[2] ).Modulus() +
1148 ( xyz[2] - xyz[0] ).Modulus() ) / 3;
1149 gp_XYZ gc = ( xyz[0] + xyz[1] + xyz[2] ) / 3;
1150 RestrictLocalSize( ngMesh, gc, size, /*overrideMinH=*/false );
1152 ngMesh.AddSurfaceElement (tri);
1153 #ifdef DUMP_TRIANGLES
1154 cout << tri << endl;
1157 if ( isInternalFace )
1159 swap( tri[1], tri[2] );
1160 ngMesh.AddSurfaceElement (tri);
1161 #ifdef DUMP_TRIANGLES
1162 cout << tri << endl;
1167 if ( quadHelper ) // remember medium nodes of sub-meshes
1169 SMDS_ElemIteratorPtr faces = smDS->GetElements();
1170 while ( faces->more() )
1172 const SMDS_MeshElement* f = faces->next();
1173 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshFace*>( f )))
1179 } // case TopAbs_FACE
1181 case TopAbs_VERTEX: { // VERTEX
1182 // --------------------------
1183 // issue 0021405. Add node only if a VERTEX is shared by a not meshed EDGE,
1184 // else netgen removes a free node and nodeVector becomes invalid
1185 PShapeIteratorPtr ansIt = helper.GetAncestors( sm->GetSubShape(),
1189 while ( const TopoDS_Shape* e = ansIt->next() )
1191 SMESH_subMesh* eSub = helper.GetMesh()->GetSubMesh( *e );
1192 if (( toAdd = ( eSub->IsEmpty() && !SMESH_Algo::isDegenerated( TopoDS::Edge( *e )))))
1197 SMDS_NodeIteratorPtr nodeIt = smDS->GetNodes();
1198 if ( nodeIt->more() )
1199 ngNodeId( nodeIt->next(), ngMesh, nodeNgIdMap );
1205 } // loop on submeshes
1208 nodeVec.resize( ngMesh.GetNP() + 1 );
1209 TNode2IdMap::iterator node_NgId, nodeNgIdEnd = nodeNgIdMap.end();
1210 for ( node_NgId = nodeNgIdMap.begin(); node_NgId != nodeNgIdEnd; ++node_NgId)
1211 nodeVec[ node_NgId->second ] = node_NgId->first;
1216 //================================================================================
1218 * \brief Duplicate mesh faces on internal geom faces
1220 //================================================================================
1222 void NETGENPlugin_Mesher::FixIntFaces(const netgen::OCCGeometry& occgeom,
1223 netgen::Mesh& ngMesh,
1224 NETGENPlugin_Internals& internalShapes)
1226 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1228 // find ng indices of internal faces
1230 for ( int ngFaceID = 1; ngFaceID <= occgeom.fmap.Extent(); ++ngFaceID )
1232 int smeshID = meshDS->ShapeToIndex( occgeom.fmap( ngFaceID ));
1233 if ( internalShapes.isInternalShape( smeshID ))
1234 ngFaceIds.insert( ngFaceID );
1236 if ( !ngFaceIds.empty() )
1239 int i, nbFaces = ngMesh.GetNSE();
1240 for ( i = 1; i <= nbFaces; ++i)
1242 netgen::Element2d elem = ngMesh.SurfaceElement(i);
1243 if ( ngFaceIds.count( elem.GetIndex() ))
1245 swap( elem[1], elem[2] );
1246 ngMesh.AddSurfaceElement (elem);
1252 //================================================================================
1254 * \brief Tries to heal the mesh on a FACE. The FACE is supposed to be partially
1255 * meshed due to NETGEN failure
1256 * \param [in] occgeom - geometry
1257 * \param [in,out] ngMesh - the mesh to fix
1258 * \param [inout] faceID - ID of the FACE to fix the mesh on
1259 * \return bool - is mesh is or becomes OK
1261 //================================================================================
1263 bool NETGENPlugin_Mesher::FixFaceMesh(const netgen::OCCGeometry& occgeom,
1264 netgen::Mesh& ngMesh,
1267 // we address a case where the FACE is almost fully meshed except small holes
1268 // of usually triangular shape at FACE boundary (IPAL52861)
1270 // The case appeared to be not simple: holes only look triangular but
1271 // indeed are a self intersecting polygon. A reason of the bug was in coincident
1272 // NG points on a seam edge. But the code below is very nice, leave it for
1277 if ( occgeom.fmap.Extent() < faceID )
1279 //const TopoDS_Face& face = TopoDS::Face( occgeom.fmap( faceID ));
1281 // find free links on the FACE
1282 NCollection_Map<Link> linkMap;
1283 for ( int iF = 1; iF <= ngMesh.GetNSE(); ++iF )
1285 const netgen::Element2d& elem = ngMesh.SurfaceElement(iF);
1286 if ( faceID != elem.GetIndex() )
1288 int n0 = elem[ elem.GetNP() - 1 ];
1289 for ( int i = 0; i < elem.GetNP(); ++i )
1292 Link link( n0, n1 );
1293 if ( !linkMap.Add( link ))
1294 linkMap.Remove( link );
1298 // add/remove boundary links
1299 for ( int iSeg = 1; iSeg <= ngMesh.GetNSeg(); ++iSeg )
1301 const netgen::Segment& seg = ngMesh.LineSegment( iSeg );
1302 if ( seg.si != faceID ) // !edgeIDs.Contains( seg.edgenr ))
1304 Link link( seg[1], seg[0] ); // reverse!!!
1305 if ( !linkMap.Add( link ))
1306 linkMap.Remove( link );
1308 if ( linkMap.IsEmpty() )
1310 if ( linkMap.Extent() < 3 )
1313 // make triangles of the links
1315 netgen::Element2d tri(3);
1316 tri.SetIndex ( faceID );
1318 NCollection_Map<Link>::Iterator linkIt( linkMap );
1319 Link link1 = linkIt.Value();
1320 // look for a link connected to link1
1321 NCollection_Map<Link>::Iterator linkIt2 = linkIt;
1322 for ( linkIt2.Next(); linkIt2.More(); linkIt2.Next() )
1324 const Link& link2 = linkIt2.Value();
1325 if ( link2.IsConnected( link1 ))
1327 // look for a link connected to both link1 and link2
1328 NCollection_Map<Link>::Iterator linkIt3 = linkIt2;
1329 for ( linkIt3.Next(); linkIt3.More(); linkIt3.Next() )
1331 const Link& link3 = linkIt3.Value();
1332 if ( link3.IsConnected( link1 ) &&
1333 link3.IsConnected( link2 ) )
1338 tri[2] = ( link2.Contains( link1.n1 ) ? link2.n1 : link3.n1 );
1339 if ( tri[0] == tri[2] || tri[1] == tri[2] )
1341 ngMesh.AddSurfaceElement( tri );
1343 // prepare for the next tria search
1344 if ( linkMap.Extent() == 3 )
1346 linkMap.Remove( link3 );
1347 linkMap.Remove( link2 );
1349 linkMap.Remove( link1 );
1350 link1 = linkIt.Value();
1363 //================================================================================
1364 // define gp_XY_Subtracted pointer to function calling gp_XY::Subtracted(gp_XY)
1365 gp_XY_FunPtr(Subtracted);
1366 //gp_XY_FunPtr(Added);
1368 //================================================================================
1370 * \brief Evaluate distance between two 2d points along the surface
1372 //================================================================================
1374 double evalDist( const gp_XY& uv1,
1376 const Handle(Geom_Surface)& surf,
1377 const int stopHandler=-1)
1379 if ( stopHandler > 0 ) // continue recursion
1381 gp_XY mid = SMESH_MesherHelper::GetMiddleUV( surf, uv1, uv2 );
1382 return evalDist( uv1,mid, surf, stopHandler-1 ) + evalDist( mid,uv2, surf, stopHandler-1 );
1384 double dist3D = surf->Value( uv1.X(), uv1.Y() ).Distance( surf->Value( uv2.X(), uv2.Y() ));
1385 if ( stopHandler == 0 ) // stop recursion
1388 // start recursion if necessary
1389 double dist2D = SMESH_MesherHelper::ApplyIn2D(surf, uv1, uv2, gp_XY_Subtracted, 0).Modulus();
1390 if ( fabs( dist3D - dist2D ) < dist2D * 1e-10 )
1391 return dist3D; // equal parametrization of a planar surface
1393 return evalDist( uv1, uv2, surf, 3 ); // start recursion
1396 //================================================================================
1398 * \brief Data of vertex internal in geom face
1400 //================================================================================
1404 gp_XY uv; //!< UV in face parametric space
1405 int ngId; //!< ng id of corrsponding node
1406 gp_XY uvClose; //!< UV of closest boundary node
1407 int ngIdClose; //!< ng id of closest boundary node
1410 //================================================================================
1412 * \brief Data of vertex internal in solid
1414 //================================================================================
1418 int ngId; //!< ng id of corresponding node
1419 int ngIdClose; //!< ng id of closest 2d mesh element
1420 int ngIdCloseN; //!< ng id of closest node of the closest 2d mesh element
1423 inline double dist2(const netgen::MeshPoint& p1, const netgen::MeshPoint& p2)
1425 return gp_Pnt( NGPOINT_COORDS(p1)).SquareDistance( gp_Pnt( NGPOINT_COORDS(p2)));
1429 //================================================================================
1431 * \brief Make netgen take internal vertices in faces into account by adding
1432 * segments including internal vertices
1434 * This function works in supposition that 1D mesh is already computed in ngMesh
1436 //================================================================================
1438 void NETGENPlugin_Mesher::AddIntVerticesInFaces(const netgen::OCCGeometry& occgeom,
1439 netgen::Mesh& ngMesh,
1440 vector<const SMDS_MeshNode*>& nodeVec,
1441 NETGENPlugin_Internals& internalShapes)
1443 if ((int) nodeVec.size() < ngMesh.GetNP() )
1444 nodeVec.resize( ngMesh.GetNP(), 0 );
1446 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1447 SMESH_MesherHelper helper( internalShapes.getMesh() );
1449 const map<int,list<int> >& face2Vert = internalShapes.getFacesWithVertices();
1450 map<int,list<int> >::const_iterator f2v = face2Vert.begin();
1451 for ( ; f2v != face2Vert.end(); ++f2v )
1453 const TopoDS_Face& face = TopoDS::Face( meshDS->IndexToShape( f2v->first ));
1454 if ( face.IsNull() ) continue;
1455 int faceNgID = occgeom.fmap.FindIndex (face);
1456 if ( faceNgID < 0 ) continue;
1458 TopLoc_Location loc;
1459 Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc);
1461 helper.SetSubShape( face );
1462 helper.SetElementsOnShape( true );
1464 // Get data of internal vertices and add them to ngMesh
1466 multimap< double, TIntVData > dist2VData; // sort vertices by distance from boundary nodes
1468 int i, nbSegInit = ngMesh.GetNSeg();
1470 // boundary characteristics
1471 double totSegLen2D = 0;
1474 const list<int>& iVertices = f2v->second;
1475 list<int>::const_iterator iv = iVertices.begin();
1476 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1479 // get node on vertex
1480 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1481 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1484 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1485 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1486 nV = SMESH_Algo::VertexNode( V, meshDS );
1487 if ( !nV ) continue;
1490 netgen::MeshPoint mp( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1491 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1492 vData.ngId = ngMesh.GetNP();
1493 nodeVec.push_back( nV );
1497 vData.uv = helper.GetNodeUV( face, nV, 0, &uvOK );
1498 if ( !uvOK ) helper.CheckNodeUV( face, nV, vData.uv, BRep_Tool::Tolerance(V),/*force=*/1);
1500 // loop on all segments of the face to find the node closest to vertex and to count
1501 // average segment 2d length
1502 double closeDist2 = numeric_limits<double>::max(), dist2;
1504 for (i = 1; i <= ngMesh.GetNSeg(); ++i)
1506 netgen::Segment & seg = ngMesh.LineSegment(i);
1507 if ( seg.si != faceNgID ) continue;
1509 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1511 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1512 if ( ngIdLast == seg[ iEnd ] ) continue;
1513 dist2 = helper.ApplyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1514 if ( dist2 < closeDist2 )
1515 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1516 ngIdLast = seg[ iEnd ];
1520 totSegLen2D += helper.ApplyIn2D(surf, uv[0], uv[1], gp_XY_Subtracted, false).Modulus();
1524 dist2VData.insert( make_pair( closeDist2, vData ));
1527 if ( totNbSeg == 0 ) break;
1528 double avgSegLen2d = totSegLen2D / totNbSeg;
1530 // Loop on vertices to add segments
1532 multimap< double, TIntVData >::iterator dist_vData = dist2VData.begin();
1533 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1535 double closeDist2 = dist_vData->first, dist2;
1536 TIntVData & vData = dist_vData->second;
1538 // try to find more close node among segments added for internal vertices
1539 for (i = nbSegInit+1; i <= ngMesh.GetNSeg(); ++i)
1541 netgen::Segment & seg = ngMesh.LineSegment(i);
1542 if ( seg.si != faceNgID ) continue;
1544 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1546 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1547 dist2 = helper.ApplyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1548 if ( dist2 < closeDist2 )
1549 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1552 // decide whether to use the closest node as the second end of segment or to
1553 // create a new point
1554 int segEnd1 = vData.ngId;
1555 int segEnd2 = vData.ngIdClose; // to use closest node
1556 gp_XY uvV = vData.uv, uvP = vData.uvClose;
1557 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1558 double nodeDist2D = sqrt( closeDist2 );
1559 double nodeDist3D = evalDist( vData.uv, vData.uvClose, surf );
1560 bool avgLenOK = ( avgSegLen2d < 0.75 * nodeDist2D );
1561 bool hintLenOK = ( segLenHint < 0.75 * nodeDist3D );
1562 //cout << "uvV " << uvV.X() <<","<<uvV.Y() << " ";
1563 if ( hintLenOK || avgLenOK )
1565 // create a point between the closest node and V
1568 double r = min( 0.5, ( hintLenOK ? segLenHint/nodeDist3D : avgSegLen2d/nodeDist2D ));
1569 // direction from V to closet node in 2D
1570 gp_Dir2d v2n( helper.ApplyIn2D(surf, uvP, uvV, gp_XY_Subtracted, false ));
1572 uvP = vData.uv + r * nodeDist2D * v2n.XY();
1573 gp_Pnt P = surf->Value( uvP.X(), uvP.Y() ).Transformed( loc );
1575 netgen::MeshPoint mp( netgen::Point<3> (P.X(), P.Y(), P.Z()));
1576 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1577 segEnd2 = ngMesh.GetNP();
1578 //cout << "Middle " << r << " uv " << uvP.X() << "," << uvP.Y() << "( " << ngMesh.Point(segEnd2).X()<<","<<ngMesh.Point(segEnd2).Y()<<","<<ngMesh.Point(segEnd2).Z()<<" )"<< endl;
1579 SMDS_MeshNode * nP = helper.AddNode(P.X(), P.Y(), P.Z());
1580 nodeVec.push_back( nP );
1582 //else cout << "at Node " << " uv " << uvP.X() << "," << uvP.Y() << endl;
1585 netgen::Segment seg;
1587 if ( segEnd1 > segEnd2 ) swap( segEnd1, segEnd2 ), swap( uvV, uvP );
1588 seg[0] = segEnd1; // ng node id
1589 seg[1] = segEnd2; // ng node id
1590 seg.edgenr = ngMesh.GetNSeg() + 1;// segment id
1593 seg.epgeominfo[ 0 ].dist = 0; // param on curve
1594 seg.epgeominfo[ 0 ].u = uvV.X();
1595 seg.epgeominfo[ 0 ].v = uvV.Y();
1596 seg.epgeominfo[ 1 ].dist = 1; // param on curve
1597 seg.epgeominfo[ 1 ].u = uvP.X();
1598 seg.epgeominfo[ 1 ].v = uvP.Y();
1600 // seg.epgeominfo[ 0 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1601 // seg.epgeominfo[ 1 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1603 ngMesh.AddSegment (seg);
1605 // add reverse segment
1606 swap( seg[0], seg[1] );
1607 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1608 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1609 ngMesh.AddSegment (seg);
1615 //================================================================================
1617 * \brief Make netgen take internal vertices in solids into account by adding
1618 * faces including internal vertices
1620 * This function works in supposition that 2D mesh is already computed in ngMesh
1622 //================================================================================
1624 void NETGENPlugin_Mesher::AddIntVerticesInSolids(const netgen::OCCGeometry& occgeom,
1625 netgen::Mesh& ngMesh,
1626 vector<const SMDS_MeshNode*>& nodeVec,
1627 NETGENPlugin_Internals& internalShapes)
1629 #ifdef DUMP_TRIANGLES_SCRIPT
1630 // create a python script making a mesh containing triangles added for internal vertices
1631 ofstream py(DUMP_TRIANGLES_SCRIPT);
1632 py << "import SMESH"<< endl
1633 << "from salome.smesh import smeshBuilder"<<endl
1634 << "smesh = smeshBuilder.New(salome.myStudy)"<<endl
1635 << "m = smesh.Mesh(name='triangles')" << endl;
1637 if ((int) nodeVec.size() < ngMesh.GetNP() )
1638 nodeVec.resize( ngMesh.GetNP(), 0 );
1640 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1641 SMESH_MesherHelper helper( internalShapes.getMesh() );
1643 const map<int,list<int> >& so2Vert = internalShapes.getSolidsWithVertices();
1644 map<int,list<int> >::const_iterator s2v = so2Vert.begin();
1645 for ( ; s2v != so2Vert.end(); ++s2v )
1647 const TopoDS_Shape& solid = meshDS->IndexToShape( s2v->first );
1648 if ( solid.IsNull() ) continue;
1649 int solidNgID = occgeom.somap.FindIndex (solid);
1650 if ( solidNgID < 0 && !occgeom.somap.IsEmpty() ) continue;
1652 helper.SetSubShape( solid );
1653 helper.SetElementsOnShape( true );
1655 // find ng indices of faces within the solid
1657 for (TopExp_Explorer fExp(solid, TopAbs_FACE); fExp.More(); fExp.Next() )
1658 ngFaceIds.insert( occgeom.fmap.FindIndex( fExp.Current() ));
1659 if ( ngFaceIds.size() == 1 && *ngFaceIds.begin() == 0 )
1660 ngFaceIds.insert( 1 );
1662 // Get data of internal vertices and add them to ngMesh
1664 multimap< double, TIntVSoData > dist2VData; // sort vertices by distance from ng faces
1666 int i, nbFaceInit = ngMesh.GetNSE();
1668 // boundary characteristics
1669 double totSegLen = 0;
1672 const list<int>& iVertices = s2v->second;
1673 list<int>::const_iterator iv = iVertices.begin();
1674 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1677 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1679 // get node on vertex
1680 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1683 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1684 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1685 nV = SMESH_Algo::VertexNode( V, meshDS );
1686 if ( !nV ) continue;
1689 netgen::MeshPoint mpV( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1690 ngMesh.AddPoint ( mpV, 1, netgen::FIXEDPOINT );
1691 vData.ngId = ngMesh.GetNP();
1692 nodeVec.push_back( nV );
1694 // loop on all 2d elements to find the one closest to vertex and to count
1695 // average segment length
1696 double closeDist2 = numeric_limits<double>::max(), avgDist2;
1697 for (i = 1; i <= ngMesh.GetNSE(); ++i)
1699 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1700 if ( !ngFaceIds.count( elem.GetIndex() )) continue;
1702 multimap< double, int> dist2nID; // sort nodes of element by distance from V
1703 for ( int j = 0; j < elem.GetNP(); ++j)
1705 netgen::MeshPoint mp = ngMesh.Point( elem[j] );
1706 double d2 = dist2( mpV, mp );
1707 dist2nID.insert( make_pair( d2, elem[j] ));
1708 avgDist2 += d2 / elem.GetNP();
1710 totNbSeg++, totSegLen+= sqrt( dist2( mp, ngMesh.Point( elem[(j+1)%elem.GetNP()])));
1712 double dist = dist2nID.begin()->first; //avgDist2;
1713 if ( dist < closeDist2 )
1714 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= dist;
1716 dist2VData.insert( make_pair( closeDist2, vData ));
1719 if ( totNbSeg == 0 ) break;
1720 double avgSegLen = totSegLen / totNbSeg;
1722 // Loop on vertices to add triangles
1724 multimap< double, TIntVSoData >::iterator dist_vData = dist2VData.begin();
1725 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1727 double closeDist2 = dist_vData->first;
1728 TIntVSoData & vData = dist_vData->second;
1730 const netgen::MeshPoint& mpV = ngMesh.Point( vData.ngId );
1732 // try to find more close face among ones added for internal vertices
1733 for (i = nbFaceInit+1; i <= ngMesh.GetNSE(); ++i)
1735 double avgDist2 = 0;
1736 multimap< double, int> dist2nID;
1737 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1738 for ( int j = 0; j < elem.GetNP(); ++j)
1740 double d = dist2( mpV, ngMesh.Point( elem[j] ));
1741 dist2nID.insert( make_pair( d, elem[j] ));
1742 avgDist2 += d / elem.GetNP();
1743 if ( avgDist2 < closeDist2 )
1744 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= avgDist2;
1747 // sort nodes of the closest face by angle with vector from V to the closest node
1748 const double tol = numeric_limits<double>::min();
1749 map< double, int > angle2ID;
1750 const netgen::Element2d& closeFace = ngMesh.SurfaceElement( vData.ngIdClose );
1751 netgen::MeshPoint mp[2];
1752 mp[0] = ngMesh.Point( vData.ngIdCloseN );
1753 gp_XYZ p1( NGPOINT_COORDS( mp[0] ));
1754 gp_XYZ pV( NGPOINT_COORDS( mpV ));
1755 gp_Vec v2p1( pV, p1 );
1756 double distN1 = v2p1.Magnitude();
1757 if ( distN1 <= tol ) continue;
1759 for ( int j = 0; j < closeFace.GetNP(); ++j)
1761 mp[1] = ngMesh.Point( closeFace[j] );
1762 gp_Vec v2p( pV, gp_Pnt( NGPOINT_COORDS( mp[1] )) );
1763 angle2ID.insert( make_pair( v2p1.Angle( v2p ), closeFace[j]));
1765 // get node with angle of 60 degrees or greater
1766 map< double, int >::iterator angle_id = angle2ID.lower_bound( 60. * M_PI / 180. );
1767 if ( angle_id == angle2ID.end() ) angle_id = --angle2ID.end();
1768 const double minAngle = 30. * M_PI / 180.;
1769 const double angle = angle_id->first;
1770 bool angleOK = ( angle > minAngle );
1772 // find points to create a triangle
1773 netgen::Element2d tri(3);
1775 tri[0] = vData.ngId;
1776 tri[1] = vData.ngIdCloseN; // to use the closest nodes
1777 tri[2] = angle_id->second; // to use the node with best angle
1779 // decide whether to use the closest node and the node with best angle or to create new ones
1780 for ( int isBestAngleN = 0; isBestAngleN < 2; ++isBestAngleN )
1782 bool createNew = !angleOK; //, distOK = true;
1784 int triInd = isBestAngleN ? 2 : 1;
1785 mp[isBestAngleN] = ngMesh.Point( tri[triInd] );
1790 double distN2 = sqrt( dist2( mpV, mp[isBestAngleN]));
1791 createNew = ( fabs( distN2 - distN1 ) > 0.25 * distN1 );
1793 else if ( angle < tol )
1795 v2p1.SetX( v2p1.X() + 1e-3 );
1801 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1802 bool avgLenOK = ( avgSegLen < 0.75 * distN1 );
1803 bool hintLenOK = ( segLenHint < 0.75 * distN1 );
1804 createNew = (createNew || avgLenOK || hintLenOK );
1805 // we create a new node not closer than 0.5 to the closest face
1806 // in order not to clash with other close face
1807 double r = min( 0.5, ( hintLenOK ? segLenHint : avgSegLen ) / distN1 );
1808 distFromV = r * distN1;
1812 // create a new point, between the node and the vertex if angleOK
1813 gp_XYZ p( NGPOINT_COORDS( mp[isBestAngleN] ));
1814 gp_Vec v2p( pV, p ); v2p.Normalize();
1815 if ( isBestAngleN && !angleOK )
1816 p = p1 + gp_Dir( v2p.XYZ() - v2p1.XYZ()).XYZ() * distN1 * 0.95;
1818 p = pV + v2p.XYZ() * distFromV;
1820 if ( !isBestAngleN ) p1 = p, distN1 = distFromV;
1822 mp[isBestAngleN].SetPoint( netgen::Point<3> (p.X(), p.Y(), p.Z()));
1823 ngMesh.AddPoint ( mp[isBestAngleN], 1, netgen::SURFACEPOINT );
1824 tri[triInd] = ngMesh.GetNP();
1825 nodeVec.push_back( helper.AddNode( p.X(), p.Y(), p.Z()) );
1828 ngMesh.AddSurfaceElement (tri);
1829 swap( tri[1], tri[2] );
1830 ngMesh.AddSurfaceElement (tri);
1832 #ifdef DUMP_TRIANGLES_SCRIPT
1833 py << "n1 = m.AddNode( "<< mpV(0)<<", "<< mpV(1)<<", "<< mpV(2)<<") "<< endl
1834 << "n2 = m.AddNode( "<< mp[0](0)<<", "<< mp[0](1)<<", "<< mp[0](2)<<") "<< endl
1835 << "n3 = m.AddNode( "<< mp[1](0)<<", "<< mp[1](1)<<", "<< mp[1](2)<<" )" << endl
1836 << "m.AddFace([n1,n2,n3])" << endl;
1838 } // loop on internal vertices of a solid
1840 } // loop on solids with internal vertices
1843 //================================================================================
1845 * \brief Fill netgen mesh with segments of a FACE
1846 * \param ngMesh - netgen mesh
1847 * \param geom - container of OCCT geometry to mesh
1848 * \param wires - data of nodes on FACE boundary
1849 * \param helper - mesher helper holding the FACE
1850 * \param nodeVec - vector of nodes in which node index == netgen ID
1851 * \retval SMESH_ComputeErrorPtr - error description
1853 //================================================================================
1855 SMESH_ComputeErrorPtr
1856 NETGENPlugin_Mesher::AddSegmentsToMesh(netgen::Mesh& ngMesh,
1857 netgen::OCCGeometry& geom,
1858 const TSideVector& wires,
1859 SMESH_MesherHelper& helper,
1860 vector< const SMDS_MeshNode* > & nodeVec,
1861 const bool overrideMinH)
1863 // ----------------------------
1864 // Check wires and count nodes
1865 // ----------------------------
1867 for ( size_t iW = 0; iW < wires.size(); ++iW )
1869 StdMeshers_FaceSidePtr wire = wires[ iW ];
1870 if ( wire->MissVertexNode() )
1872 // Commented for issue 0020960. It worked for the case, let's wait for case where it doesn't.
1873 // It seems that there is no reason for this limitation
1875 // (new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH, "Missing nodes on vertices"));
1877 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1878 if ((int) uvPtVec.size() != wire->NbPoints() )
1879 return SMESH_ComputeError::New(COMPERR_BAD_INPUT_MESH,
1880 SMESH_Comment("Unexpected nb of points on wire ") << iW
1881 << ": " << uvPtVec.size()<<" != "<<wire->NbPoints());
1882 nbNodes += wire->NbPoints();
1884 nodeVec.reserve( nodeVec.size() + nbNodes + 1 );
1885 if ( nodeVec.empty() )
1886 nodeVec.push_back( 0 );
1888 // -----------------
1890 // -----------------
1892 const bool wasNgMeshEmpty = ( ngMesh.GetNP() < 1 ); /* true => this method is called by
1893 NETGENPlugin_NETGEN_2D_ONLY */
1895 // map for nodes on vertices since they can be shared between wires
1896 // ( issue 0020676, face_int_box.brep) and nodes built by NETGEN
1897 map<const SMDS_MeshNode*, int > node2ngID;
1898 if ( !wasNgMeshEmpty ) // fill node2ngID with nodes built by NETGEN
1900 set< int > subIDs; // ids of sub-shapes of the FACE
1901 for ( size_t iW = 0; iW < wires.size(); ++iW )
1903 StdMeshers_FaceSidePtr wire = wires[ iW ];
1904 for ( int iE = 0, nbE = wire->NbEdges(); iE < nbE; ++iE )
1906 subIDs.insert( wire->EdgeID( iE ));
1907 subIDs.insert( helper.GetMeshDS()->ShapeToIndex( wire->FirstVertex( iE )));
1910 for ( size_t ngID = 1; ngID < nodeVec.size(); ++ngID )
1911 if ( subIDs.count( nodeVec[ngID]->getshapeId() ))
1912 node2ngID.insert( make_pair( nodeVec[ngID], ngID ));
1915 const int solidID = 0, faceID = geom.fmap.FindIndex( helper.GetSubShape() );
1916 if ( ngMesh.GetNFD() < 1 )
1917 ngMesh.AddFaceDescriptor( netgen::FaceDescriptor( faceID, solidID, solidID, 0 ));
1919 for ( size_t iW = 0; iW < wires.size(); ++iW )
1921 StdMeshers_FaceSidePtr wire = wires[ iW ];
1922 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1923 const int nbSegments = wire->NbPoints() - 1;
1925 // assure the 1st node to be in node2ngID, which is needed to correctly
1926 // "close chain of segments" (see below) in case if the 1st node is not
1927 // onVertex because it is on a Viscous layer
1928 node2ngID.insert( make_pair( uvPtVec[ 0 ].node, ngMesh.GetNP() + 1 ));
1930 // compute length of every segment
1931 vector<double> segLen( nbSegments );
1932 for ( int i = 0; i < nbSegments; ++i )
1933 segLen[i] = SMESH_TNodeXYZ( uvPtVec[ i ].node ).Distance( uvPtVec[ i+1 ].node );
1935 int edgeID = 1, posID = -2;
1936 bool isInternalWire = false;
1937 double vertexNormPar = 0;
1938 const int prevNbNGSeg = ngMesh.GetNSeg();
1939 for ( int i = 0; i < nbSegments; ++i ) // loop on segments
1941 // Add the first point of a segment
1943 const SMDS_MeshNode * n = uvPtVec[ i ].node;
1944 const int posShapeID = n->getshapeId();
1945 bool onVertex = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX );
1946 bool onEdge = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_EDGE );
1948 // skip nodes on degenerated edges
1949 if ( helper.IsDegenShape( posShapeID ) &&
1950 helper.IsDegenShape( uvPtVec[ i+1 ].node->getshapeId() ))
1953 int ngID1 = ngMesh.GetNP() + 1, ngID2 = ngID1+1;
1954 if ( onVertex || ( !wasNgMeshEmpty && onEdge ) || helper.IsRealSeam( posShapeID ))
1955 ngID1 = node2ngID.insert( make_pair( n, ngID1 )).first->second;
1956 if ( ngID1 > ngMesh.GetNP() )
1958 netgen::MeshPoint mp( netgen::Point<3> (n->X(), n->Y(), n->Z()) );
1959 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1960 nodeVec.push_back( n );
1962 else // n is in ngMesh already, and ngID2 in prev segment is wrong
1964 ngID2 = ngMesh.GetNP() + 1;
1965 if ( i > 0 ) // prev segment belongs to same wire
1967 netgen::Segment& prevSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
1974 netgen::Segment seg;
1976 seg[0] = ngID1; // ng node id
1977 seg[1] = ngID2; // ng node id
1978 seg.edgenr = ngMesh.GetNSeg() + 1; // ng segment id
1979 seg.si = faceID; // = geom.fmap.FindIndex (face);
1981 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1983 const UVPtStruct& pnt = uvPtVec[ i + iEnd ];
1985 seg.epgeominfo[ iEnd ].dist = pnt.param; // param on curve
1986 seg.epgeominfo[ iEnd ].u = pnt.u;
1987 seg.epgeominfo[ iEnd ].v = pnt.v;
1989 // find out edge id and node parameter on edge
1990 onVertex = ( pnt.normParam + 1e-10 > vertexNormPar );
1991 if ( onVertex || posShapeID != posID )
1994 double normParam = pnt.normParam;
1996 normParam = 0.5 * ( uvPtVec[ i ].normParam + uvPtVec[ i+1 ].normParam );
1997 int edgeIndexInWire = wire->EdgeIndex( normParam );
1998 vertexNormPar = wire->LastParameter( edgeIndexInWire );
1999 const TopoDS_Edge& edge = wire->Edge( edgeIndexInWire );
2000 edgeID = geom.emap.FindIndex( edge );
2002 isInternalWire = ( edge.Orientation() == TopAbs_INTERNAL );
2003 // if ( onVertex ) // param on curve is different on each of two edges
2004 // seg.epgeominfo[ iEnd ].dist = helper.GetNodeU( edge, pnt.node );
2006 seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
2009 ngMesh.AddSegment (seg);
2011 // restrict size of elements near the segment
2012 SMESH_TNodeXYZ np1( n ), np2( uvPtVec[ i+1 ].node );
2013 // get an average size of adjacent segments to avoid sharp change of
2014 // element size (regression on issue 0020452, note 0010898)
2015 int iPrev = SMESH_MesherHelper::WrapIndex( i-1, nbSegments );
2016 int iNext = SMESH_MesherHelper::WrapIndex( i+1, nbSegments );
2017 double sumH = segLen[ iPrev ] + segLen[ i ] + segLen[ iNext ];
2018 int nbSeg = ( int( segLen[ iPrev ] > sumH / 100.) +
2019 int( segLen[ i ] > sumH / 100.) +
2020 int( segLen[ iNext ] > sumH / 100.));
2022 RestrictLocalSize( ngMesh, 0.5*(np1+np2), sumH / nbSeg, overrideMinH );
2024 if ( isInternalWire )
2026 swap (seg[0], seg[1]);
2027 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
2028 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
2029 ngMesh.AddSegment (seg);
2031 } // loop on segments on a wire
2033 // close chain of segments
2034 if ( nbSegments > 0 )
2036 netgen::Segment& lastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() - int( isInternalWire ));
2037 const SMDS_MeshNode * lastNode = uvPtVec.back().node;
2038 lastSeg[1] = node2ngID.insert( make_pair( lastNode, lastSeg[1] )).first->second;
2039 if ( lastSeg[1] > ngMesh.GetNP() )
2041 netgen::MeshPoint mp( netgen::Point<3> (lastNode->X(), lastNode->Y(), lastNode->Z()) );
2042 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
2043 nodeVec.push_back( lastNode );
2045 if ( isInternalWire )
2047 netgen::Segment& realLastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
2048 realLastSeg[0] = lastSeg[1];
2052 #ifdef DUMP_SEGMENTS
2053 cout << "BEGIN WIRE " << iW << endl;
2054 for ( int i = prevNbNGSeg+1; i <= ngMesh.GetNSeg(); ++i )
2056 netgen::Segment& seg = ngMesh.LineSegment( i );
2058 netgen::Segment& prevSeg = ngMesh.LineSegment( i-1 );
2059 if ( seg[0] == prevSeg[1] && seg[1] == prevSeg[0] )
2061 cout << "Segment: " << seg.edgenr << endl << "\tis REVRESE of the previous one" << endl;
2065 cout << "Segment: " << seg.edgenr << endl
2066 << "\tp1: " << seg[0] << " n" << nodeVec[ seg[0]]->GetID() << endl
2067 << "\tp2: " << seg[1] << " n" << nodeVec[ seg[1]]->GetID() << endl
2068 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
2069 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
2070 << "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
2071 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
2072 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl
2073 << "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
2075 cout << "--END WIRE " << iW << endl;
2077 SMESH_Comment __not_unused_variable( prevNbNGSeg );
2080 } // loop on WIREs of a FACE
2082 // add a segment instead of an internal vertex
2083 if ( wasNgMeshEmpty )
2085 NETGENPlugin_Internals intShapes( *helper.GetMesh(), helper.GetSubShape(), /*is3D=*/false );
2086 AddIntVerticesInFaces( geom, ngMesh, nodeVec, intShapes );
2088 ngMesh.CalcSurfacesOfNode();
2093 //================================================================================
2095 * \brief Fill SMESH mesh according to contents of netgen mesh
2096 * \param occgeo - container of OCCT geometry to mesh
2097 * \param ngMesh - netgen mesh
2098 * \param initState - bn of entities in netgen mesh before computing
2099 * \param sMesh - SMESH mesh to fill in
2100 * \param nodeVec - vector of nodes in which node index == netgen ID
2101 * \param comment - returns problem description
2102 * \param quadHelper - holder of medium nodes of sub-meshes
2103 * \retval int - error
2105 //================================================================================
2107 int NETGENPlugin_Mesher::FillSMesh(const netgen::OCCGeometry& occgeo,
2108 netgen::Mesh& ngMesh,
2109 const NETGENPlugin_ngMeshInfo& initState,
2111 std::vector<const SMDS_MeshNode*>& nodeVec,
2112 SMESH_Comment& comment,
2113 SMESH_MesherHelper* quadHelper)
2115 int nbNod = ngMesh.GetNP();
2116 int nbSeg = ngMesh.GetNSeg();
2117 int nbFac = ngMesh.GetNSE();
2118 int nbVol = ngMesh.GetNE();
2120 SMESHDS_Mesh* meshDS = sMesh.GetMeshDS();
2122 // quadHelper is used for either
2123 // 1) making quadratic elements when a lower dimention mesh is loaded
2124 // to SMESH before convertion to quadratic by NETGEN
2125 // 2) sewing of quadratic elements with quadratic elements of sub-meshes
2126 if ( quadHelper && !quadHelper->GetIsQuadratic() && quadHelper->GetTLinkNodeMap().empty() )
2129 // -------------------------------------
2130 // Create and insert nodes into nodeVec
2131 // -------------------------------------
2133 nodeVec.resize( nbNod + 1 );
2134 int i, nbInitNod = initState._nbNodes;
2135 for (i = nbInitNod+1; i <= nbNod; ++i )
2137 const netgen::MeshPoint& ngPoint = ngMesh.Point(i);
2138 SMDS_MeshNode* node = NULL;
2139 TopoDS_Vertex aVert;
2140 // First, netgen creates nodes on vertices in occgeo.vmap,
2141 // so node index corresponds to vertex index
2142 // but (issue 0020776) netgen does not create nodes with equal coordinates
2143 if ( i-nbInitNod <= occgeo.vmap.Extent() )
2145 gp_Pnt p ( NGPOINT_COORDS(ngPoint) );
2146 for (int iV = i-nbInitNod; aVert.IsNull() && iV <= occgeo.vmap.Extent(); ++iV)
2148 aVert = TopoDS::Vertex( occgeo.vmap( iV ));
2149 gp_Pnt pV = BRep_Tool::Pnt( aVert );
2150 if ( p.SquareDistance( pV ) > 1e-20 )
2153 node = const_cast<SMDS_MeshNode*>( SMESH_Algo::VertexNode( aVert, meshDS ));
2156 if (!node) // node not found on vertex
2158 node = meshDS->AddNode( NGPOINT_COORDS( ngPoint ));
2159 if (!aVert.IsNull())
2160 meshDS->SetNodeOnVertex(node, aVert);
2165 // -------------------------------------------
2166 // Create mesh segments along geometric edges
2167 // -------------------------------------------
2169 int nbInitSeg = initState._nbSegments;
2170 for (i = nbInitSeg+1; i <= nbSeg; ++i )
2172 const netgen::Segment& seg = ngMesh.LineSegment(i);
2174 int pinds[3] = { seg.pnums[0], seg.pnums[1], seg.pnums[2] };
2177 for (int j=0; j < 3; ++j)
2179 int pind = pinds[j];
2180 if (pind <= 0 || !nodeVec_ACCESS(pind))
2188 int aGeomEdgeInd = seg.epgeominfo[j].edgenr;
2189 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
2190 aEdge = TopoDS::Edge(occgeo.emap(aGeomEdgeInd));
2192 param = seg.epgeominfo[j].dist;
2195 else // middle point
2197 param = param2 * 0.5;
2199 if (!aEdge.IsNull() && nodeVec_ACCESS(pind)->getshapeId() < 1)
2201 meshDS->SetNodeOnEdge(nodeVec_ACCESS(pind), aEdge, param);
2206 SMDS_MeshEdge* edge = 0;
2207 if (nbp == 2) // second order ?
2209 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1])))
2211 if ( quadHelper ) // final mesh must be quadratic
2212 edge = quadHelper->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
2214 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
2218 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
2219 nodeVec_ACCESS(pinds[2])))
2221 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
2222 nodeVec_ACCESS(pinds[2]));
2226 if ( comment.empty() ) comment << "Cannot create a mesh edge";
2227 MESSAGE("Cannot create a mesh edge");
2228 nbSeg = nbFac = nbVol = 0;
2231 if ( !aEdge.IsNull() && edge->getshapeId() < 1 )
2232 meshDS->SetMeshElementOnShape(edge, aEdge);
2234 else if ( comment.empty() )
2236 comment << "Invalid netgen segment #" << i;
2240 // ----------------------------------------
2241 // Create mesh faces along geometric faces
2242 // ----------------------------------------
2244 int nbInitFac = initState._nbFaces;
2245 int quadFaceID = ngMesh.GetNFD() + 1;
2246 if ( nbInitFac < nbFac )
2247 // add a faces descriptor to exclude qudrangle elements generated by NETGEN
2248 // from computation of 3D mesh
2249 ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(quadFaceID, /*solid1=*/0, /*solid2=*/0, 0));
2251 vector<const SMDS_MeshNode*> nodes;
2252 for (i = nbInitFac+1; i <= nbFac; ++i )
2254 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
2255 const int aGeomFaceInd = elem.GetIndex();
2257 if (aGeomFaceInd > 0 && aGeomFaceInd <= occgeo.fmap.Extent())
2258 aFace = TopoDS::Face(occgeo.fmap(aGeomFaceInd));
2260 for ( int j = 1; j <= elem.GetNP(); ++j )
2262 int pind = elem.PNum(j);
2263 if ( pind < 1 || pind >= (int) nodeVec.size() )
2265 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind))
2267 nodes.push_back( node );
2268 if (!aFace.IsNull() && node->getshapeId() < 1)
2270 const netgen::PointGeomInfo& pgi = elem.GeomInfoPi(j);
2271 meshDS->SetNodeOnFace(node, aFace, pgi.u, pgi.v);
2275 if ((int) nodes.size() != elem.GetNP() )
2277 if ( comment.empty() )
2278 comment << "Invalid netgen 2d element #" << i;
2279 continue; // bad node ids
2281 SMDS_MeshFace* face = NULL;
2282 switch (elem.GetType())
2285 if ( quadHelper ) // final mesh must be quadratic
2286 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2]);
2288 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2]);
2291 if ( quadHelper ) // final mesh must be quadratic
2292 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2294 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2295 // exclude qudrangle elements from computation of 3D mesh
2296 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2299 nodes[5] = mediumNode( nodes[0],nodes[1],nodes[5], quadHelper );
2300 nodes[3] = mediumNode( nodes[1],nodes[2],nodes[3], quadHelper );
2301 nodes[4] = mediumNode( nodes[2],nodes[0],nodes[4], quadHelper );
2302 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[5],nodes[3],nodes[4]);
2305 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2306 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2307 nodes[5] = mediumNode( nodes[2],nodes[3],nodes[5], quadHelper );
2308 nodes[6] = mediumNode( nodes[3],nodes[0],nodes[6], quadHelper );
2309 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3],
2310 nodes[4],nodes[7],nodes[5],nodes[6]);
2311 // exclude qudrangle elements from computation of 3D mesh
2312 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2315 MESSAGE("NETGEN created a face of unexpected type, ignoring");
2320 if ( comment.empty() ) comment << "Cannot create a mesh face";
2321 MESSAGE("Cannot create a mesh face");
2322 nbSeg = nbFac = nbVol = 0;
2325 if ( !aFace.IsNull() )
2326 meshDS->SetMeshElementOnShape( face, aFace );
2329 // ------------------
2330 // Create tetrahedra
2331 // ------------------
2333 for ( i = 1; i <= nbVol; ++i )
2335 const netgen::Element& elem = ngMesh.VolumeElement(i);
2336 int aSolidInd = elem.GetIndex();
2337 TopoDS_Solid aSolid;
2338 if ( aSolidInd > 0 && aSolidInd <= occgeo.somap.Extent() )
2339 aSolid = TopoDS::Solid(occgeo.somap(aSolidInd));
2341 for ( int j = 1; j <= elem.GetNP(); ++j )
2343 int pind = elem.PNum(j);
2344 if ( pind < 1 || pind >= (int)nodeVec.size() )
2346 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind) )
2348 nodes.push_back(node);
2349 if ( !aSolid.IsNull() && node->getshapeId() < 1 )
2350 meshDS->SetNodeInVolume(node, aSolid);
2353 if ((int) nodes.size() != elem.GetNP() )
2355 if ( comment.empty() )
2356 comment << "Invalid netgen 3d element #" << i;
2359 SMDS_MeshVolume* vol = NULL;
2360 switch ( elem.GetType() )
2363 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3]);
2366 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2367 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2368 nodes[5] = mediumNode( nodes[2],nodes[0],nodes[5], quadHelper );
2369 nodes[6] = mediumNode( nodes[0],nodes[3],nodes[6], quadHelper );
2370 nodes[8] = mediumNode( nodes[1],nodes[3],nodes[8], quadHelper );
2371 nodes[9] = mediumNode( nodes[2],nodes[3],nodes[9], quadHelper );
2372 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3],
2373 nodes[4],nodes[7],nodes[5],nodes[6],nodes[8],nodes[9]);
2376 MESSAGE("NETGEN created a volume of unexpected type, ignoring");
2381 if ( comment.empty() ) comment << "Cannot create a mesh volume";
2382 MESSAGE("Cannot create a mesh volume");
2383 nbSeg = nbFac = nbVol = 0;
2386 if (!aSolid.IsNull())
2387 meshDS->SetMeshElementOnShape(vol, aSolid);
2389 return comment.empty() ? 0 : 1;
2394 //================================================================================
2396 * \brief Convert error into text
2398 //================================================================================
2400 std::string text(int err)
2405 SMESH_Comment("Error in netgen::OCCGenerateMesh() at ") << netgen::multithread.task;
2408 //================================================================================
2410 * \brief Convert exception into text
2412 //================================================================================
2414 std::string text(Standard_Failure& ex)
2416 SMESH_Comment str("Exception in netgen::OCCGenerateMesh()");
2417 str << " at " << netgen::multithread.task
2418 << ": " << ex.DynamicType()->Name();
2419 if ( ex.GetMessageString() && strlen( ex.GetMessageString() ))
2420 str << ": " << ex.GetMessageString();
2423 //================================================================================
2425 * \brief Convert exception into text
2427 //================================================================================
2429 std::string text(netgen::NgException& ex)
2431 SMESH_Comment str("NgException");
2432 if ( strlen( netgen::multithread.task ) > 0 )
2433 str << " at " << netgen::multithread.task;
2434 str << ": " << ex.What();
2438 //================================================================================
2440 * \brief Looks for triangles lying on a SOLID
2442 //================================================================================
2444 bool hasBadElemOnSolid( const list<const SMDS_MeshElement*>& elems,
2445 SMESH_subMesh* solidSM )
2447 TopTools_IndexedMapOfShape solidSubs;
2448 TopExp::MapShapes( solidSM->GetSubShape(), solidSubs );
2449 SMESHDS_Mesh* mesh = solidSM->GetFather()->GetMeshDS();
2451 list<const SMDS_MeshElement*>::const_iterator e = elems.begin();
2452 for ( ; e != elems.end(); ++e )
2454 const SMDS_MeshElement* elem = *e;
2455 // if ( elem->GetType() != SMDSAbs_Face ) -- 23047
2457 int nbNodesOnSolid = 0, nbNodes = elem->NbNodes();
2458 SMDS_NodeIteratorPtr nIt = elem->nodeIterator();
2459 while ( nIt->more() )
2461 const SMDS_MeshNode* n = nIt->next();
2462 const TopoDS_Shape& s = mesh->IndexToShape( n->getshapeId() );
2463 nbNodesOnSolid += ( !s.IsNull() && solidSubs.Contains( s ));
2464 if ( nbNodesOnSolid > 2 ||
2465 nbNodesOnSolid == nbNodes)
2472 const double edgeMeshingTime = 0.001;
2473 const double faceMeshingTime = 0.019;
2474 const double edgeFaceMeshingTime = edgeMeshingTime + faceMeshingTime;
2475 const double faceOptimizTime = 0.06;
2476 const double voluMeshingTime = 0.15;
2477 const double volOptimizeTime = 0.77;
2480 //=============================================================================
2482 * Here we are going to use the NETGEN mesher
2484 //=============================================================================
2486 bool NETGENPlugin_Mesher::Compute()
2488 NETGENPlugin_NetgenLibWrapper ngLib;
2490 netgen::MeshingParameters& mparams = netgen::mparam;
2492 SMESH_ComputeErrorPtr error = SMESH_ComputeError::New();
2493 SMESH_MesherHelper quadHelper( *_mesh );
2494 quadHelper.SetIsQuadratic( mparams.secondorder );
2496 static string debugFile = "/tmp/ngMesh.py"; /* to call toPython( _ngMesh, debugFile )
2497 while debugging netgen */
2498 // -------------------------
2499 // Prepare OCC geometry
2500 // -------------------------
2502 netgen::OCCGeometry occgeo;
2503 list< SMESH_subMesh* > meshedSM[3]; // for 0-2 dimensions
2504 NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
2505 PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
2508 _totalTime = edgeFaceMeshingTime;
2510 _totalTime += faceOptimizTime;
2512 _totalTime += voluMeshingTime + ( _optimize ? volOptimizeTime : 0 );
2513 double doneTime = 0;
2516 _curShapeIndex = -1;
2518 // -------------------------
2519 // Generate the mesh
2520 // -------------------------
2523 NETGENPlugin_ngMeshInfo initState; // it remembers size of ng mesh equal to size of Smesh
2525 SMESH_Comment comment;
2528 // vector of nodes in which node index == netgen ID
2529 vector< const SMDS_MeshNode* > nodeVec;
2537 // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
2538 mparams.uselocalh = false;
2539 mparams.grading = 0.8; // not limitited size growth
2541 if ( _simpleHyp->GetNumberOfSegments() )
2543 mparams.maxh = occgeo.boundingbox.Diam();
2546 mparams.maxh = _simpleHyp->GetLocalLength();
2549 if ( mparams.maxh == 0.0 )
2550 mparams.maxh = occgeo.boundingbox.Diam();
2551 if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
2552 mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
2554 // Local size on faces
2555 occgeo.face_maxh = mparams.maxh;
2557 // Let netgen create _ngMesh and calculate element size on not meshed shapes
2561 int startWith = netgen::MESHCONST_ANALYSE;
2562 int endWith = netgen::MESHCONST_ANALYSE;
2567 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2569 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2571 if(netgen::multithread.terminate)
2574 comment << text(err);
2576 catch (Standard_Failure& ex)
2578 comment << text(ex);
2580 catch (netgen::NgException & ex)
2582 comment << text(ex);
2583 if ( mparams.meshsizefilename )
2584 throw SMESH_ComputeError(COMPERR_BAD_PARMETERS, comment );
2586 err = 0; //- MESHCONST_ANALYSE isn't so important step
2589 ngLib.setMesh(( Ng_Mesh*) _ngMesh );
2591 _ngMesh->ClearFaceDescriptors(); // we make descriptors our-self
2593 if ( !mparams.uselocalh ) // mparams.grading is not taken into account yet
2594 _ngMesh->LocalHFunction().SetGrading( mparams.grading );
2598 // Pass 1D simple parameters to NETGEN
2599 // --------------------------------
2600 int nbSeg = _simpleHyp->GetNumberOfSegments();
2601 double segSize = _simpleHyp->GetLocalLength();
2602 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
2604 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
2606 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
2607 setLocalSize( e, segSize, *_ngMesh );
2610 else // if ( ! _simpleHyp )
2612 // Local size on shapes
2613 SetLocalSize( occgeo, *_ngMesh );
2616 // Precompute internal edges (issue 0020676) in order to
2617 // add mesh on them correctly (twice) to netgen mesh
2618 if ( !err && internals.hasInternalEdges() )
2620 // load internal shapes into OCCGeometry
2621 netgen::OCCGeometry intOccgeo;
2622 internals.getInternalEdges( intOccgeo.fmap, intOccgeo.emap, intOccgeo.vmap, meshedSM );
2623 intOccgeo.boundingbox = occgeo.boundingbox;
2624 intOccgeo.shape = occgeo.shape;
2625 intOccgeo.face_maxh.SetSize(intOccgeo.fmap.Extent());
2626 intOccgeo.face_maxh = netgen::mparam.maxh;
2627 netgen::Mesh *tmpNgMesh = NULL;
2631 // compute local H on internal shapes in the main mesh
2632 //OCCSetLocalMeshSize(intOccgeo, *_ngMesh); it deletes _ngMesh->localH
2634 // let netgen create a temporary mesh
2636 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2638 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2640 if(netgen::multithread.terminate)
2643 // copy LocalH from the main to temporary mesh
2644 initState.transferLocalH( _ngMesh, tmpNgMesh );
2646 // compute mesh on internal edges
2647 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2649 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2651 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2653 comment << text(err);
2655 catch (Standard_Failure& ex)
2657 comment << text(ex);
2660 initState.restoreLocalH( tmpNgMesh );
2662 // fill SMESH by netgen mesh
2663 vector< const SMDS_MeshNode* > tmpNodeVec;
2664 FillSMesh( intOccgeo, *tmpNgMesh, initState, *_mesh, tmpNodeVec, comment );
2665 err = ( err || !comment.empty() );
2667 nglib::Ng_DeleteMesh((nglib::Ng_Mesh*)tmpNgMesh);
2670 // Fill _ngMesh with nodes and segments of computed submeshes
2673 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_0D ]) &&
2674 FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_1D ], &quadHelper));
2676 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2681 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2686 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2688 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2690 if(netgen::multithread.terminate)
2693 comment << text(err);
2695 catch (Standard_Failure& ex)
2697 comment << text(ex);
2702 _ticTime = ( doneTime += edgeMeshingTime ) / _totalTime / _progressTic;
2704 mparams.uselocalh = true; // restore as it is used at surface optimization
2706 // ---------------------
2707 // compute surface mesh
2708 // ---------------------
2711 // Pass 2D simple parameters to NETGEN
2713 if ( double area = _simpleHyp->GetMaxElementArea() ) {
2715 mparams.maxh = sqrt(2. * area/sqrt(3.0));
2716 mparams.grading = 0.4; // moderate size growth
2719 // length from edges
2720 if ( _ngMesh->GetNSeg() ) {
2721 double edgeLength = 0;
2722 TopTools_MapOfShape visitedEdges;
2723 for ( TopExp_Explorer exp( _shape, TopAbs_EDGE ); exp.More(); exp.Next() )
2724 if( visitedEdges.Add(exp.Current()) )
2725 edgeLength += SMESH_Algo::EdgeLength( TopoDS::Edge( exp.Current() ));
2726 // we have to multiply length by 2 since for each TopoDS_Edge there
2727 // are double set of NETGEN edges, in other words, we have to
2728 // divide _ngMesh->GetNSeg() by 2.
2729 mparams.maxh = 2*edgeLength / _ngMesh->GetNSeg();
2732 mparams.maxh = 1000;
2734 mparams.grading = 0.2; // slow size growth
2736 mparams.quad = _simpleHyp->GetAllowQuadrangles();
2737 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2738 _ngMesh->SetGlobalH (mparams.maxh);
2739 netgen::Box<3> bb = occgeo.GetBoundingBox();
2740 bb.Increase (bb.Diam()/20);
2741 _ngMesh->SetLocalH (bb.PMin(), bb.PMax(), mparams.grading);
2744 // Care of vertices internal in faces (issue 0020676)
2745 if ( internals.hasInternalVertexInFace() )
2747 // store computed segments in SMESH in order not to create SMESH
2748 // edges for ng segments added by AddIntVerticesInFaces()
2749 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2750 // add segments to faces with internal vertices
2751 AddIntVerticesInFaces( occgeo, *_ngMesh, nodeVec, internals );
2752 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2755 // Build viscous layers
2756 if ( _isViscousLayers2D ||
2757 StdMeshers_ViscousLayers2D::HasProxyMesh( TopoDS::Face( occgeo.fmap(1) ), *_mesh ))
2759 if ( !internals.hasInternalVertexInFace() ) {
2760 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2761 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2763 SMESH_ProxyMesh::Ptr viscousMesh;
2764 SMESH_MesherHelper helper( *_mesh );
2765 for ( int faceID = 1; faceID <= occgeo.fmap.Extent(); ++faceID )
2767 const TopoDS_Face& F = TopoDS::Face( occgeo.fmap( faceID ));
2768 viscousMesh = StdMeshers_ViscousLayers2D::Compute( *_mesh, F );
2771 if ( viscousMesh->NbProxySubMeshes() == 0 )
2773 // exclude from computation ng segments built on EDGEs of F
2774 for (int i = 1; i <= _ngMesh->GetNSeg(); i++)
2776 netgen::Segment & seg = _ngMesh->LineSegment(i);
2777 if (seg.si == faceID)
2780 // add new segments to _ngMesh instead of excluded ones
2781 helper.SetSubShape( F );
2783 StdMeshers_FaceSide::GetFaceWires( F, *_mesh, /*skipMediumNodes=*/true,
2784 error, viscousMesh );
2785 error = AddSegmentsToMesh( *_ngMesh, occgeo, wires, helper, nodeVec );
2787 if ( !error ) error = SMESH_ComputeError::New();
2789 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2792 // Let netgen compute 2D mesh
2793 startWith = netgen::MESHCONST_MESHSURFACE;
2794 endWith = _optimize ? netgen::MESHCONST_OPTSURFACE : netgen::MESHCONST_MESHSURFACE;
2799 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2801 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2803 if(netgen::multithread.terminate)
2806 comment << text (err);
2808 catch (Standard_Failure& ex)
2810 comment << text(ex);
2811 //err = 1; -- try to make volumes anyway
2813 catch (netgen::NgException exc)
2815 comment << text(exc);
2816 //err = 1; -- try to make volumes anyway
2821 doneTime += faceMeshingTime + ( _optimize ? faceOptimizTime : 0 );
2822 _ticTime = doneTime / _totalTime / _progressTic;
2824 // ---------------------
2825 // generate volume mesh
2826 // ---------------------
2827 // Fill _ngMesh with nodes and faces of computed 2D submeshes
2828 if ( !err && _isVolume && ( !meshedSM[ MeshDim_2D ].empty() || mparams.quad ))
2830 // load SMESH with computed segments and faces
2831 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2833 // compute pyramids on quadrangles
2834 SMESH_ProxyMesh::Ptr proxyMesh;
2835 if ( _mesh->NbQuadrangles() > 0 )
2836 for ( int iS = 1; iS <= occgeo.somap.Extent(); ++iS )
2838 StdMeshers_QuadToTriaAdaptor* Adaptor = new StdMeshers_QuadToTriaAdaptor;
2839 proxyMesh.reset( Adaptor );
2841 int nbPyrams = _mesh->NbPyramids();
2842 Adaptor->Compute( *_mesh, occgeo.somap(iS) );
2843 if ( nbPyrams != _mesh->NbPyramids() )
2845 list< SMESH_subMesh* > quadFaceSM;
2846 for (TopExp_Explorer face(occgeo.somap(iS), TopAbs_FACE); face.More(); face.Next())
2847 if ( Adaptor->GetProxySubMesh( face.Current() ))
2849 quadFaceSM.push_back( _mesh->GetSubMesh( face.Current() ));
2850 meshedSM[ MeshDim_2D ].remove( quadFaceSM.back() );
2852 FillNgMesh(occgeo, *_ngMesh, nodeVec, quadFaceSM, &quadHelper, proxyMesh);
2855 // fill _ngMesh with faces of sub-meshes
2856 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_2D ], &quadHelper));
2857 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2858 //toPython( _ngMesh, "/tmp/ngPython.py");
2860 if (!err && _isVolume)
2862 // Pass 3D simple parameters to NETGEN
2863 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
2864 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
2866 if ( double vol = simple3d->GetMaxElementVolume() ) {
2868 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
2869 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2872 // length from faces
2873 mparams.maxh = _ngMesh->AverageH();
2875 _ngMesh->SetGlobalH (mparams.maxh);
2876 mparams.grading = 0.4;
2878 _ngMesh->CalcLocalH(mparams.grading);
2880 _ngMesh->CalcLocalH();
2883 // Care of vertices internal in solids and internal faces (issue 0020676)
2884 if ( internals.hasInternalVertexInSolid() || internals.hasInternalFaces() )
2886 // store computed faces in SMESH in order not to create SMESH
2887 // faces for ng faces added here
2888 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2889 // add ng faces to solids with internal vertices
2890 AddIntVerticesInSolids( occgeo, *_ngMesh, nodeVec, internals );
2891 // duplicate mesh faces on internal faces
2892 FixIntFaces( occgeo, *_ngMesh, internals );
2893 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2895 // Let netgen compute 3D mesh
2896 startWith = endWith = netgen::MESHCONST_MESHVOLUME;
2901 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2903 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2905 if(netgen::multithread.terminate)
2908 if ( comment.empty() ) // do not overwrite a previos error
2909 comment << text(err);
2911 catch (Standard_Failure& ex)
2913 if ( comment.empty() ) // do not overwrite a previos error
2914 comment << text(ex);
2917 catch (netgen::NgException exc)
2919 if ( comment.empty() ) // do not overwrite a previos error
2920 comment << text(exc);
2923 _ticTime = ( doneTime += voluMeshingTime ) / _totalTime / _progressTic;
2925 // Let netgen optimize 3D mesh
2926 if ( !err && _optimize )
2928 startWith = endWith = netgen::MESHCONST_OPTVOLUME;
2933 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2935 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2937 if(netgen::multithread.terminate)
2940 if ( comment.empty() ) // do not overwrite a previos error
2941 comment << text(err);
2943 catch (Standard_Failure& ex)
2945 if ( comment.empty() ) // do not overwrite a previos error
2946 comment << text(ex);
2948 catch (netgen::NgException exc)
2950 if ( comment.empty() ) // do not overwrite a previos error
2951 comment << text(exc);
2955 if (!err && mparams.secondorder > 0)
2960 if ( !meshedSM[ MeshDim_1D ].empty() )
2962 // remove segments not attached to geometry (IPAL0052479)
2963 for (int i = 1; i <= _ngMesh->GetNSeg(); ++i)
2965 const netgen::Segment & seg = _ngMesh->LineSegment (i);
2966 if ( seg.epgeominfo[ 0 ].edgenr == 0 )
2967 _ngMesh->DeleteSegment( i );
2969 _ngMesh->Compress();
2971 // convert to quadratic
2972 netgen::OCCRefinementSurfaces ref (occgeo);
2973 ref.MakeSecondOrder (*_ngMesh);
2975 // care of elements already loaded to SMESH
2976 // if ( initState._nbSegments > 0 )
2977 // makeQuadratic( occgeo.emap, _mesh );
2978 // if ( initState._nbFaces > 0 )
2979 // makeQuadratic( occgeo.fmap, _mesh );
2981 catch (Standard_Failure& ex)
2983 if ( comment.empty() ) // do not overwrite a previos error
2984 comment << "Exception in netgen at passing to 2nd order ";
2986 catch (netgen::NgException exc)
2988 if ( comment.empty() ) // do not overwrite a previos error
2989 comment << exc.What();
2994 _ticTime = 0.98 / _progressTic;
2996 //int nbNod = _ngMesh->GetNP();
2997 //int nbSeg = _ngMesh->GetNSeg();
2998 int nbFac = _ngMesh->GetNSE();
2999 int nbVol = _ngMesh->GetNE();
3000 bool isOK = ( !err && (_isVolume ? (nbVol > 0) : (nbFac > 0)) );
3002 // Feed back the SMESHDS with the generated Nodes and Elements
3003 if ( true /*isOK*/ ) // get whatever built
3005 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
3007 if ( quadHelper.GetIsQuadratic() ) // remove free nodes
3008 for ( size_t i = 0; i < nodeVec.size(); ++i )
3009 if ( nodeVec[i] && nodeVec[i]->NbInverseElements() == 0 )
3010 _mesh->GetMeshDS()->RemoveFreeNode( nodeVec[i], 0, /*fromGroups=*/false );
3012 SMESH_ComputeErrorPtr readErr = ReadErrors(nodeVec);
3013 if ( readErr && !readErr->myBadElements.empty() )
3016 if ( !comment.empty() && !readErr->myComment.empty() ) comment += "\n";
3017 comment += readErr->myComment;
3019 if ( error->IsOK() && ( !isOK || comment.size() > 0 ))
3020 error->myName = COMPERR_ALGO_FAILED;
3021 if ( !comment.empty() )
3022 error->myComment = comment;
3024 // SetIsAlwaysComputed( true ) to empty sub-meshes, which
3025 // appear if the geometry contains coincident sub-shape due
3026 // to bool merge_solids = 1; in netgen/libsrc/occ/occgenmesh.cpp
3027 const int nbMaps = 2;
3028 const TopTools_IndexedMapOfShape* geoMaps[nbMaps] =
3029 { & occgeo.vmap, & occgeo.emap/*, & occgeo.fmap*/ };
3030 for ( int iMap = 0; iMap < nbMaps; ++iMap )
3031 for (int i = 1; i <= geoMaps[iMap]->Extent(); i++)
3032 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( geoMaps[iMap]->FindKey(i)))
3033 if ( !sm->IsMeshComputed() )
3034 sm->SetIsAlwaysComputed( true );
3036 // set bad compute error to subshapes of all failed sub-shapes
3037 if ( !error->IsOK() )
3039 bool pb2D = false, pb3D = false;
3040 for (int i = 1; i <= occgeo.fmap.Extent(); i++) {
3041 int status = occgeo.facemeshstatus[i-1];
3042 if (status == netgen::FACE_MESHED_OK ) continue;
3043 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.fmap( i ))) {
3044 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
3045 if ( !smError || smError->IsOK() ) {
3046 if ( status == netgen::FACE_FAILED )
3047 smError.reset( new SMESH_ComputeError( *error ));
3049 smError.reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED, "Ignored" ));
3050 if ( SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
3051 smError->myName = COMPERR_WARNING;
3053 pb2D = pb2D || smError->IsKO();
3056 if ( !pb2D ) // all faces are OK
3057 for (int i = 1; i <= occgeo.somap.Extent(); i++)
3058 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.somap( i )))
3060 bool smComputed = nbVol && !sm->IsEmpty();
3061 if ( smComputed && internals.hasInternalVertexInSolid( sm->GetId() ))
3063 int nbIntV = internals.getSolidsWithVertices().find( sm->GetId() )->second.size();
3064 SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
3065 smComputed = ( smDS->NbElements() > 0 || smDS->NbNodes() > nbIntV );
3067 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
3068 if ( !smComputed && ( !smError || smError->IsOK() ))
3070 smError.reset( new SMESH_ComputeError( *error ));
3071 if ( nbVol && SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
3073 smError->myName = COMPERR_WARNING;
3075 else if ( !smError->myBadElements.empty() ) // bad surface mesh
3077 if ( !hasBadElemOnSolid( smError->myBadElements, sm ))
3081 pb3D = pb3D || ( smError && smError->IsKO() );
3083 if ( !pb2D && !pb3D )
3084 err = 0; // no fatal errors, only warnings
3087 ngLib._isComputeOk = !err;
3092 //=============================================================================
3096 //=============================================================================
3097 bool NETGENPlugin_Mesher::Evaluate(MapShapeNbElems& aResMap)
3099 netgen::MeshingParameters& mparams = netgen::mparam;
3102 // -------------------------
3103 // Prepare OCC geometry
3104 // -------------------------
3105 netgen::OCCGeometry occgeo;
3106 list< SMESH_subMesh* > meshedSM[4]; // for 0-3 dimensions
3107 NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
3108 PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
3110 bool tooManyElems = false;
3111 const int hugeNb = std::numeric_limits<int>::max() / 100;
3116 // pass 1D simple parameters to NETGEN
3119 // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
3120 mparams.uselocalh = false;
3121 mparams.grading = 0.8; // not limitited size growth
3123 if ( _simpleHyp->GetNumberOfSegments() )
3125 mparams.maxh = occgeo.boundingbox.Diam();
3128 mparams.maxh = _simpleHyp->GetLocalLength();
3131 if ( mparams.maxh == 0.0 )
3132 mparams.maxh = occgeo.boundingbox.Diam();
3133 if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
3134 mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
3136 // let netgen create _ngMesh and calculate element size on not meshed shapes
3137 NETGENPlugin_NetgenLibWrapper ngLib;
3138 netgen::Mesh *ngMesh = NULL;
3142 int startWith = netgen::MESHCONST_ANALYSE;
3143 int endWith = netgen::MESHCONST_MESHEDGES;
3145 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, mparams, startWith, endWith);
3147 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, startWith, endWith, optstr);
3150 if(netgen::multithread.terminate)
3153 ngLib.setMesh(( Ng_Mesh*) ngMesh );
3155 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( _shape ))
3156 sm->GetComputeError().reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED ));
3161 // Pass 1D simple parameters to NETGEN
3162 // --------------------------------
3163 int nbSeg = _simpleHyp->GetNumberOfSegments();
3164 double segSize = _simpleHyp->GetLocalLength();
3165 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
3167 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
3169 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
3170 setLocalSize( e, segSize, *ngMesh );
3173 else // if ( ! _simpleHyp )
3175 // Local size on shapes
3176 SetLocalSize( occgeo, *ngMesh );
3178 // calculate total nb of segments and length of edges
3179 double fullLen = 0.0;
3181 int entity = mparams.secondorder > 0 ? SMDSEntity_Quad_Edge : SMDSEntity_Edge;
3182 TopTools_DataMapOfShapeInteger Edge2NbSeg;
3183 for (TopExp_Explorer exp(_shape, TopAbs_EDGE); exp.More(); exp.Next())
3185 TopoDS_Edge E = TopoDS::Edge( exp.Current() );
3186 if( !Edge2NbSeg.Bind(E,0) )
3189 double aLen = SMESH_Algo::EdgeLength(E);
3192 vector<int>& aVec = aResMap[_mesh->GetSubMesh(E)];
3194 aVec.resize( SMDSEntity_Last, 0);
3196 fullNbSeg += aVec[ entity ];
3199 // store nb of segments computed by Netgen
3200 NCollection_Map<Link> linkMap;
3201 for (int i = 1; i <= ngMesh->GetNSeg(); ++i )
3203 const netgen::Segment& seg = ngMesh->LineSegment(i);
3204 Link link(seg[0], seg[1]);
3205 if ( !linkMap.Add( link )) continue;
3206 int aGeomEdgeInd = seg.epgeominfo[0].edgenr;
3207 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
3209 vector<int>& aVec = aResMap[_mesh->GetSubMesh(occgeo.emap(aGeomEdgeInd))];
3213 // store nb of nodes on edges computed by Netgen
3214 TopTools_DataMapIteratorOfDataMapOfShapeInteger Edge2NbSegIt(Edge2NbSeg);
3215 for (; Edge2NbSegIt.More(); Edge2NbSegIt.Next())
3217 vector<int>& aVec = aResMap[_mesh->GetSubMesh(Edge2NbSegIt.Key())];
3218 if ( aVec[ entity ] > 1 && aVec[ SMDSEntity_Node ] == 0 )
3219 aVec[SMDSEntity_Node] = mparams.secondorder > 0 ? 2*aVec[ entity ]-1 : aVec[ entity ]-1;
3221 fullNbSeg += aVec[ entity ];
3222 Edge2NbSeg( Edge2NbSegIt.Key() ) = aVec[ entity ];
3224 if ( fullNbSeg == 0 )
3231 if ( double area = _simpleHyp->GetMaxElementArea() ) {
3233 mparams.maxh = sqrt(2. * area/sqrt(3.0));
3234 mparams.grading = 0.4; // moderate size growth
3237 // length from edges
3238 mparams.maxh = fullLen/fullNbSeg;
3239 mparams.grading = 0.2; // slow size growth
3242 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3243 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3245 for (TopExp_Explorer exp(_shape, TopAbs_FACE); exp.More(); exp.Next())
3247 TopoDS_Face F = TopoDS::Face( exp.Current() );
3248 SMESH_subMesh *sm = _mesh->GetSubMesh(F);
3250 BRepGProp::SurfaceProperties(F,G);
3251 double anArea = G.Mass();
3252 tooManyElems = tooManyElems || ( anArea/hugeNb > mparams.maxh*mparams.maxh );
3254 if ( !tooManyElems )
3256 TopTools_MapOfShape egdes;
3257 for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next())
3258 if ( egdes.Add( exp1.Current() ))
3259 nb1d += Edge2NbSeg.Find(exp1.Current());
3261 int nbFaces = tooManyElems ? hugeNb : int( 4*anArea / (mparams.maxh*mparams.maxh*sqrt(3.)));
3262 int nbNodes = tooManyElems ? hugeNb : (( nbFaces*3 - (nb1d-1)*2 ) / 6 + 1 );
3264 vector<int> aVec(SMDSEntity_Last, 0);
3265 if( mparams.secondorder > 0 ) {
3266 int nb1d_in = (nbFaces*3 - nb1d) / 2;
3267 aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
3268 aVec[SMDSEntity_Quad_Triangle] = nbFaces;
3271 aVec[SMDSEntity_Node] = Max ( nbNodes, 0 );
3272 aVec[SMDSEntity_Triangle] = nbFaces;
3274 aResMap[sm].swap(aVec);
3281 // pass 3D simple parameters to NETGEN
3282 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
3283 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
3285 if ( double vol = simple3d->GetMaxElementVolume() ) {
3287 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
3288 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3291 // using previous length from faces
3293 mparams.grading = 0.4;
3294 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3297 BRepGProp::VolumeProperties(_shape,G);
3298 double aVolume = G.Mass();
3299 double tetrVol = 0.1179*mparams.maxh*mparams.maxh*mparams.maxh;
3300 tooManyElems = tooManyElems || ( aVolume/hugeNb > tetrVol );
3301 int nbVols = tooManyElems ? hugeNb : int(aVolume/tetrVol);
3302 int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
3303 vector<int> aVec(SMDSEntity_Last, 0 );
3304 if ( tooManyElems ) // avoid FPE
3306 aVec[SMDSEntity_Node] = hugeNb;
3307 aVec[ mparams.secondorder > 0 ? SMDSEntity_Quad_Tetra : SMDSEntity_Tetra] = hugeNb;
3311 if( mparams.secondorder > 0 ) {
3312 aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
3313 aVec[SMDSEntity_Quad_Tetra] = nbVols;
3316 aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
3317 aVec[SMDSEntity_Tetra] = nbVols;
3320 SMESH_subMesh *sm = _mesh->GetSubMesh(_shape);
3321 aResMap[sm].swap(aVec);
3327 double NETGENPlugin_Mesher::GetProgress(const SMESH_Algo* holder,
3328 const int * algoProgressTic,
3329 const double * algoProgress) const
3331 ((int&) _progressTic ) = *algoProgressTic + 1;
3333 if ( !_occgeom ) return 0;
3335 double progress = -1;
3338 if ( _ticTime < 0 && netgen::multithread.task[0] == 'O'/*Optimizing surface*/ )
3340 ((double&) _ticTime ) = edgeFaceMeshingTime / _totalTime / _progressTic;
3342 else if ( !_optimize /*&& _occgeom->fmap.Extent() > 1*/ )
3344 int doneShapeIndex = -1;
3345 while ( doneShapeIndex+1 < _occgeom->facemeshstatus.Size() &&
3346 _occgeom->facemeshstatus[ doneShapeIndex+1 ])
3348 if ( doneShapeIndex+1 != _curShapeIndex )
3350 ((int&) _curShapeIndex) = doneShapeIndex+1;
3351 double doneShapeRate = _curShapeIndex / double( _occgeom->fmap.Extent() );
3352 double doneTime = edgeMeshingTime + doneShapeRate * faceMeshingTime;
3353 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3354 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3355 // << " " << doneTime / _totalTime / _progressTic << endl;
3359 else if ( !_optimize && _occgeom->somap.Extent() > 1 )
3361 int curShapeIndex = _curShapeIndex;
3362 if ( _ngMesh->GetNE() > 0 )
3364 netgen::Element el = (*_ngMesh)[netgen::ElementIndex( _ngMesh->GetNE()-1 )];
3365 curShapeIndex = el.GetIndex();
3367 if ( curShapeIndex != _curShapeIndex )
3369 ((int&) _curShapeIndex) = curShapeIndex;
3370 double doneShapeRate = _curShapeIndex / double( _occgeom->somap.Extent() );
3371 double doneTime = edgeFaceMeshingTime + doneShapeRate * voluMeshingTime;
3372 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3373 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3374 // << " " << doneTime / _totalTime / _progressTic << endl;
3379 progress = Max( *algoProgressTic * _ticTime, *algoProgress );
3384 netgen::multithread.task[0] == 'D'/*elaunay meshing*/ &&
3385 progress > voluMeshingTime )
3387 progress = voluMeshingTime;
3388 ((double&) _ticTime) = voluMeshingTime / _totalTime / _progressTic;
3390 ((int&) *algoProgressTic )++;
3391 ((double&) *algoProgress) = progress;
3393 //cout << progress << " " << *algoProgressTic << " " << netgen::multithread.task << " "<< _ticTime << endl;
3395 return Min( progress, 0.99 );
3398 //================================================================================
3400 * \brief Read mesh entities preventing successful computation from "test.out" file
3402 //================================================================================
3404 SMESH_ComputeErrorPtr
3405 NETGENPlugin_Mesher::ReadErrors(const vector<const SMDS_MeshNode* >& nodeVec)
3407 SMESH_ComputeErrorPtr err = SMESH_ComputeError::New
3408 (COMPERR_BAD_INPUT_MESH, "Some edges multiple times in surface mesh");
3409 SMESH_File file("test.out");
3411 vector<int> three1(3), three2(3);
3412 const char* badEdgeStr = " multiple times in surface mesh";
3413 const int badEdgeStrLen = strlen( badEdgeStr );
3414 const int nbNodes = nodeVec.size();
3416 while( !file.eof() )
3418 if ( strncmp( file, "Edge ", 5 ) == 0 &&
3419 file.getInts( two ) &&
3420 strncmp( file, badEdgeStr, badEdgeStrLen ) == 0 &&
3421 two[0] < nbNodes && two[1] < nbNodes )
3423 err->myBadElements.push_back( new SMDS_LinearEdge( nodeVec[ two[0]], nodeVec[ two[1]] ));
3424 file += badEdgeStrLen;
3426 else if ( strncmp( file, "Intersecting: ", 14 ) == 0 )
3429 // openelement 18 with open element 126
3433 const char* pos = file;
3434 bool ok = ( strncmp( file, "openelement ", 12 ) == 0 );
3435 ok = ok && file.getInts( two );
3436 ok = ok && file.getInts( three1 );
3437 ok = ok && file.getInts( three2 );
3438 for ( int i = 0; ok && i < 3; ++i )
3439 ok = ( three1[i] < nbNodes && nodeVec[ three1[i]]);
3440 for ( int i = 0; ok && i < 3; ++i )
3441 ok = ( three2[i] < nbNodes && nodeVec[ three2[i]]);
3444 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three1[0]],
3445 nodeVec[ three1[1]],
3446 nodeVec[ three1[2]]));
3447 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three2[0]],
3448 nodeVec[ three2[1]],
3449 nodeVec[ three2[2]]));
3450 err->myComment = "Intersecting triangles";
3464 size_t nbBadElems = err->myBadElements.size();
3465 if ( nbBadElems ) nbBadElems++; // avoid warning: variable set but not used
3471 //================================================================================
3473 * \brief Write a python script creating an equivalent SALOME mesh.
3474 * This is useful to see what mesh is passed as input for the next step of mesh
3475 * generation (of mesh of higher dimension)
3477 //================================================================================
3479 void NETGENPlugin_Mesher::toPython( const netgen::Mesh* ngMesh )
3481 const char* pyFile = "/tmp/ngMesh.py";
3482 ofstream outfile( pyFile, ios::out );
3483 if ( !outfile ) return;
3485 outfile << "import SMESH" << endl
3486 << "from salome.smesh import smeshBuilder" << endl
3487 << "smesh = smeshBuilder.New(salome.myStudy)" << endl
3488 << "mesh = smesh.Mesh()" << endl << endl;
3490 using namespace netgen;
3492 for (pi = PointIndex::BASE;
3493 pi < ngMesh->GetNP()+PointIndex::BASE; pi++)
3495 outfile << "mesh.AddNode( ";
3496 outfile << (*ngMesh)[pi](0) << ", ";
3497 outfile << (*ngMesh)[pi](1) << ", ";
3498 outfile << (*ngMesh)[pi](2) << ") ## "<< pi << endl;
3501 int nbDom = ngMesh->GetNDomains();
3502 for ( int i = 0; i < nbDom; ++i )
3503 outfile<< "grp" << i+1 << " = mesh.CreateEmptyGroup( SMESH.FACE, 'domain"<< i+1 << "')"<< endl;
3505 SurfaceElementIndex sei;
3506 for (sei = 0; sei < ngMesh->GetNSE(); sei++)
3508 outfile << "mesh.AddFace([ ";
3509 Element2d sel = (*ngMesh)[sei];
3510 for (int j = 0; j < sel.GetNP(); j++)
3511 outfile << sel[j] << ( j+1 < sel.GetNP() ? ", " : " ])");
3512 if ( sel.IsDeleted() ) outfile << " ## IsDeleted ";
3515 if ((*ngMesh)[sei].GetIndex())
3517 if ( int dom1 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainIn())
3518 outfile << "grp"<< dom1 <<".Add([ " << (int)sei+1 << " ])" << endl;
3519 if ( int dom2 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainOut())
3520 outfile << "grp"<< dom2 <<".Add([ " << (int)sei+1 << " ])" << endl;
3524 for (ElementIndex ei = 0; ei < ngMesh->GetNE(); ei++)
3526 Element el = (*ngMesh)[ei];
3527 outfile << "mesh.AddVolume([ ";
3528 for (int j = 0; j < el.GetNP(); j++)
3529 outfile << el[j] << ( j+1 < el.GetNP() ? ", " : " ])");
3533 for (int i = 1; i <= ngMesh->GetNSeg(); i++)
3535 const Segment & seg = ngMesh->LineSegment (i);
3536 outfile << "mesh.AddEdge([ "
3538 << seg[1] << " ])" << endl;
3540 cout << "Write " << pyFile << endl;
3543 //================================================================================
3545 * \brief Constructor of NETGENPlugin_ngMeshInfo
3547 //================================================================================
3549 NETGENPlugin_ngMeshInfo::NETGENPlugin_ngMeshInfo( netgen::Mesh* ngMesh):
3554 _nbNodes = ngMesh->GetNP();
3555 _nbSegments = ngMesh->GetNSeg();
3556 _nbFaces = ngMesh->GetNSE();
3557 _nbVolumes = ngMesh->GetNE();
3561 _nbNodes = _nbSegments = _nbFaces = _nbVolumes = 0;
3565 //================================================================================
3567 * \brief Copy LocalH member from one netgen mesh to another
3569 //================================================================================
3571 void NETGENPlugin_ngMeshInfo::transferLocalH( netgen::Mesh* fromMesh,
3572 netgen::Mesh* toMesh )
3574 if ( !fromMesh->LocalHFunctionGenerated() ) return;
3575 if ( !toMesh->LocalHFunctionGenerated() )
3577 toMesh->CalcLocalH(netgen::mparam.grading);
3579 toMesh->CalcLocalH();
3582 const size_t size = sizeof( netgen::LocalH );
3583 _copyOfLocalH = new char[ size ];
3584 memcpy( (void*)_copyOfLocalH, (void*)&toMesh->LocalHFunction(), size );
3585 memcpy( (void*)&toMesh->LocalHFunction(), (void*)&fromMesh->LocalHFunction(), size );
3588 //================================================================================
3590 * \brief Restore LocalH member of a netgen mesh
3592 //================================================================================
3594 void NETGENPlugin_ngMeshInfo::restoreLocalH( netgen::Mesh* toMesh )
3596 if ( _copyOfLocalH )
3598 const size_t size = sizeof( netgen::LocalH );
3599 memcpy( (void*)&toMesh->LocalHFunction(), (void*)_copyOfLocalH, size );
3600 delete [] _copyOfLocalH;
3605 //================================================================================
3607 * \brief Find "internal" sub-shapes
3609 //================================================================================
3611 NETGENPlugin_Internals::NETGENPlugin_Internals( SMESH_Mesh& mesh,
3612 const TopoDS_Shape& shape,
3614 : _mesh( mesh ), _is3D( is3D )
3616 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3618 TopExp_Explorer f,e;
3619 for ( f.Init( shape, TopAbs_FACE ); f.More(); f.Next() )
3621 int faceID = meshDS->ShapeToIndex( f.Current() );
3623 // find not computed internal edges
3625 for ( e.Init( f.Current().Oriented(TopAbs_FORWARD), TopAbs_EDGE ); e.More(); e.Next() )
3626 if ( e.Current().Orientation() == TopAbs_INTERNAL )
3628 SMESH_subMesh* eSM = mesh.GetSubMesh( e.Current() );
3629 if ( eSM->IsEmpty() )
3631 _e2face.insert( make_pair( eSM->GetId(), faceID ));
3632 for ( TopoDS_Iterator v(e.Current()); v.More(); v.Next() )
3633 _e2face.insert( make_pair( meshDS->ShapeToIndex( v.Value() ), faceID ));
3637 // find internal vertices in a face
3638 set<int> intVV; // issue 0020850 where same vertex is twice in a face
3639 for ( TopoDS_Iterator fSub( f.Current() ); fSub.More(); fSub.Next())
3640 if ( fSub.Value().ShapeType() == TopAbs_VERTEX )
3642 int vID = meshDS->ShapeToIndex( fSub.Value() );
3643 if ( intVV.insert( vID ).second )
3644 _f2v[ faceID ].push_back( vID );
3649 // find internal faces and their subshapes where nodes are to be doubled
3650 // to make a crack with non-sewed borders
3652 if ( f.Current().Orientation() == TopAbs_INTERNAL )
3654 _intShapes.insert( meshDS->ShapeToIndex( f.Current() ));
3657 list< TopoDS_Shape > edges;
3658 for ( e.Init( f.Current(), TopAbs_EDGE ); e.More(); e.Next())
3659 if ( SMESH_MesherHelper::NbAncestors( e.Current(), mesh, TopAbs_FACE ) > 1 )
3661 _intShapes.insert( meshDS->ShapeToIndex( e.Current() ));
3662 edges.push_back( e.Current() );
3663 // find border faces
3664 PShapeIteratorPtr fIt =
3665 SMESH_MesherHelper::GetAncestors( edges.back(),mesh,TopAbs_FACE );
3666 while ( const TopoDS_Shape* pFace = fIt->next() )
3667 if ( !pFace->IsSame( f.Current() ))
3668 _borderFaces.insert( meshDS->ShapeToIndex( *pFace ));
3671 // we consider vertex internal if it is shared by more than one internal edge
3672 list< TopoDS_Shape >::iterator edge = edges.begin();
3673 for ( ; edge != edges.end(); ++edge )
3674 for ( TopoDS_Iterator v( *edge ); v.More(); v.Next() )
3676 set<int> internalEdges;
3677 PShapeIteratorPtr eIt =
3678 SMESH_MesherHelper::GetAncestors( v.Value(),mesh,TopAbs_EDGE );
3679 while ( const TopoDS_Shape* pEdge = eIt->next() )
3681 int edgeID = meshDS->ShapeToIndex( *pEdge );
3682 if ( isInternalShape( edgeID ))
3683 internalEdges.insert( edgeID );
3685 if ( internalEdges.size() > 1 )
3686 _intShapes.insert( meshDS->ShapeToIndex( v.Value() ));
3690 } // loop on geom faces
3692 // find vertices internal in solids
3695 for ( TopExp_Explorer so(shape, TopAbs_SOLID); so.More(); so.Next())
3697 int soID = meshDS->ShapeToIndex( so.Current() );
3698 for ( TopoDS_Iterator soSub( so.Current() ); soSub.More(); soSub.Next())
3699 if ( soSub.Value().ShapeType() == TopAbs_VERTEX )
3700 _s2v[ soID ].push_back( meshDS->ShapeToIndex( soSub.Value() ));
3705 //================================================================================
3707 * \brief Find mesh faces on non-internal geom faces sharing internal edge
3708 * some nodes of which are to be doubled to make the second border of the "crack"
3710 //================================================================================
3712 void NETGENPlugin_Internals::findBorderElements( TIDSortedElemSet & borderElems )
3714 if ( _intShapes.empty() ) return;
3716 SMESH_Mesh& mesh = const_cast<SMESH_Mesh&>(_mesh);
3717 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3719 // loop on internal geom edges
3720 set<int>::const_iterator intShapeId = _intShapes.begin();
3721 for ( ; intShapeId != _intShapes.end(); ++intShapeId )
3723 const TopoDS_Shape& s = meshDS->IndexToShape( *intShapeId );
3724 if ( s.ShapeType() != TopAbs_EDGE ) continue;
3726 // get internal and non-internal geom faces sharing the internal edge <s>
3728 set<int>::iterator bordFace = _borderFaces.end();
3729 PShapeIteratorPtr faces = SMESH_MesherHelper::GetAncestors( s, _mesh, TopAbs_FACE );
3730 while ( const TopoDS_Shape* pFace = faces->next() )
3732 int faceID = meshDS->ShapeToIndex( *pFace );
3733 if ( isInternalShape( faceID ))
3736 bordFace = _borderFaces.insert( faceID ).first;
3738 if ( bordFace == _borderFaces.end() || !intFace ) continue;
3740 // get all links of mesh faces on internal geom face sharing nodes on edge <s>
3741 set< SMESH_OrientedLink > links; //!< links of faces on internal geom face
3742 list<const SMDS_MeshElement*> suspectFaces[2]; //!< mesh faces on border geom faces
3743 int nbSuspectFaces = 0;
3744 SMESHDS_SubMesh* intFaceSM = meshDS->MeshElements( intFace );
3745 if ( !intFaceSM || intFaceSM->NbElements() == 0 ) continue;
3746 SMESH_subMeshIteratorPtr smIt = mesh.GetSubMesh( s )->getDependsOnIterator(true,true);
3747 while ( smIt->more() )
3749 SMESHDS_SubMesh* sm = smIt->next()->GetSubMeshDS();
3750 if ( !sm ) continue;
3751 SMDS_NodeIteratorPtr nIt = sm->GetNodes();
3752 while ( nIt->more() )
3754 const SMDS_MeshNode* nOnEdge = nIt->next();
3755 SMDS_ElemIteratorPtr fIt = nOnEdge->GetInverseElementIterator(SMDSAbs_Face);
3756 while ( fIt->more() )
3758 const SMDS_MeshElement* f = fIt->next();
3759 const int nbNodes = f->NbCornerNodes();
3760 if ( intFaceSM->Contains( f ))
3762 for ( int i = 0; i < nbNodes; ++i )
3763 links.insert( SMESH_OrientedLink( f->GetNode(i), f->GetNode((i+1)%nbNodes)));
3768 for ( int i = 0; i < nbNodes; ++i )
3769 nbDblNodes += isInternalShape( f->GetNode(i)->getshapeId() );
3771 suspectFaces[ nbDblNodes < 2 ].push_back( f );
3777 // suspectFaces[0] having link with same orientation as mesh faces on
3778 // the internal geom face are <borderElems>. suspectFaces[1] have
3779 // only one node on edge <s>, we decide on them later (at the 2nd loop)
3780 // by links of <borderElems> found at the 1st and 2nd loops
3781 set< SMESH_OrientedLink > borderLinks;
3782 for ( int isPostponed = 0; isPostponed < 2; ++isPostponed )
3784 list<const SMDS_MeshElement*>::iterator fIt = suspectFaces[isPostponed].begin();
3785 for ( int nbF = 0; fIt != suspectFaces[isPostponed].end(); ++fIt, ++nbF )
3787 const SMDS_MeshElement* f = *fIt;
3788 bool isBorder = false, linkFound = false, borderLinkFound = false;
3789 list< SMESH_OrientedLink > faceLinks;
3790 int nbNodes = f->NbCornerNodes();
3791 for ( int i = 0; i < nbNodes; ++i )
3793 SMESH_OrientedLink link( f->GetNode(i), f->GetNode((i+1)%nbNodes));
3794 faceLinks.push_back( link );
3797 set< SMESH_OrientedLink >::iterator foundLink = links.find( link );
3798 if ( foundLink != links.end() )
3801 isBorder = ( foundLink->_reversed == link._reversed );
3802 if ( !isBorder && !isPostponed ) break;
3803 faceLinks.pop_back();
3805 else if ( isPostponed && !borderLinkFound )
3807 foundLink = borderLinks.find( link );
3808 if ( foundLink != borderLinks.end() )
3810 borderLinkFound = true;
3811 isBorder = ( foundLink->_reversed != link._reversed );
3818 borderElems.insert( f );
3819 borderLinks.insert( faceLinks.begin(), faceLinks.end() );
3821 else if ( !linkFound && !borderLinkFound )
3823 suspectFaces[1].push_back( f );
3824 if ( nbF > 2 * nbSuspectFaces )
3825 break; // dead loop protection
3832 //================================================================================
3834 * \brief put internal shapes in maps and fill in submeshes to precompute
3836 //================================================================================
3838 void NETGENPlugin_Internals::getInternalEdges( TopTools_IndexedMapOfShape& fmap,
3839 TopTools_IndexedMapOfShape& emap,
3840 TopTools_IndexedMapOfShape& vmap,
3841 list< SMESH_subMesh* > smToPrecompute[])
3843 if ( !hasInternalEdges() ) return;
3844 map<int,int>::const_iterator ev_face = _e2face.begin();
3845 for ( ; ev_face != _e2face.end(); ++ev_face )
3847 const TopoDS_Shape& ev = _mesh.GetMeshDS()->IndexToShape( ev_face->first );
3848 const TopoDS_Shape& face = _mesh.GetMeshDS()->IndexToShape( ev_face->second );
3850 ( ev.ShapeType() == TopAbs_EDGE ? emap : vmap ).Add( ev );
3852 //cout<<"INTERNAL EDGE or VERTEX "<<ev_face->first<<" on face "<<ev_face->second<<endl;
3854 smToPrecompute[ MeshDim_1D ].push_back( _mesh.GetSubMeshContaining( ev_face->first ));
3858 //================================================================================
3860 * \brief return shapes and submeshes to be meshed and already meshed boundary submeshes
3862 //================================================================================
3864 void NETGENPlugin_Internals::getInternalFaces( TopTools_IndexedMapOfShape& fmap,
3865 TopTools_IndexedMapOfShape& emap,
3866 list< SMESH_subMesh* >& intFaceSM,
3867 list< SMESH_subMesh* >& boundarySM)
3869 if ( !hasInternalFaces() ) return;
3871 // <fmap> and <emap> are for not yet meshed shapes
3872 // <intFaceSM> is for submeshes of faces
3873 // <boundarySM> is for meshed edges and vertices
3878 set<int> shapeIDs ( _intShapes );
3879 if ( !_borderFaces.empty() )
3880 shapeIDs.insert( _borderFaces.begin(), _borderFaces.end() );
3882 set<int>::const_iterator intS = shapeIDs.begin();
3883 for ( ; intS != shapeIDs.end(); ++intS )
3885 SMESH_subMesh* sm = _mesh.GetSubMeshContaining( *intS );
3887 if ( sm->GetSubShape().ShapeType() != TopAbs_FACE ) continue;
3889 intFaceSM.push_back( sm );
3891 // add submeshes of not computed internal faces
3892 if ( !sm->IsEmpty() ) continue;
3894 SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(true,true);
3895 while ( smIt->more() )
3898 const TopoDS_Shape& s = sm->GetSubShape();
3900 if ( sm->IsEmpty() )
3903 switch ( s.ShapeType() ) {
3904 case TopAbs_FACE: fmap.Add ( s ); break;
3905 case TopAbs_EDGE: emap.Add ( s ); break;
3911 if ( s.ShapeType() != TopAbs_FACE )
3912 boundarySM.push_back( sm );
3918 //================================================================================
3920 * \brief Return true if given shape is to be precomputed in order to be correctly
3921 * added to netgen mesh
3923 //================================================================================
3925 bool NETGENPlugin_Internals::isShapeToPrecompute(const TopoDS_Shape& s)
3927 int shapeID = _mesh.GetMeshDS()->ShapeToIndex( s );
3928 switch ( s.ShapeType() ) {
3929 case TopAbs_FACE : break; //return isInternalShape( shapeID ) || isBorderFace( shapeID );
3930 case TopAbs_EDGE : return isInternalEdge( shapeID );
3931 case TopAbs_VERTEX: break;
3937 //================================================================================
3939 * \brief Return SMESH
3941 //================================================================================
3943 SMESH_Mesh& NETGENPlugin_Internals::getMesh() const
3945 return const_cast<SMESH_Mesh&>( _mesh );
3948 //================================================================================
3950 * \brief Access to a counter of NETGENPlugin_NetgenLibWrapper instances
3952 //================================================================================
3954 int& NETGENPlugin_NetgenLibWrapper::instanceCounter()
3956 static int theCouner = 0;
3960 //================================================================================
3962 * \brief Initialize netgen library
3964 //================================================================================
3966 NETGENPlugin_NetgenLibWrapper::NETGENPlugin_NetgenLibWrapper()
3968 if ( instanceCounter() == 0 )
3971 ++instanceCounter();
3973 _isComputeOk = false;
3977 if ( !getenv( "KEEP_NETGEN_OUTPUT" ))
3979 // redirect all netgen output (mycout,myerr,cout) to _outputFileName
3980 _outputFileName = getOutputFileName();
3981 _ngcout = netgen::mycout;
3982 _ngcerr = netgen::myerr;
3983 netgen::mycout = new ofstream ( _outputFileName.c_str() );
3984 netgen::myerr = netgen::mycout;
3985 _coutBuffer = std::cout.rdbuf();
3987 cout << "NOTE: netgen output is redirected to file " << _outputFileName << endl;
3989 std::cout.rdbuf( netgen::mycout->rdbuf() );
3993 _ngMesh = Ng_NewMesh();
3996 //================================================================================
3998 * \brief Finish using netgen library
4000 //================================================================================
4002 NETGENPlugin_NetgenLibWrapper::~NETGENPlugin_NetgenLibWrapper()
4004 --instanceCounter();
4006 Ng_DeleteMesh( _ngMesh );
4010 std::cout.rdbuf( _coutBuffer );
4017 //================================================================================
4019 * \brief Set netgen mesh to delete at destruction
4021 //================================================================================
4023 void NETGENPlugin_NetgenLibWrapper::setMesh( Ng_Mesh* mesh )
4026 Ng_DeleteMesh( _ngMesh );
4030 //================================================================================
4032 * \brief Return a unique file name
4034 //================================================================================
4036 std::string NETGENPlugin_NetgenLibWrapper::getOutputFileName()
4038 std::string aTmpDir = SALOMEDS_Tool::GetTmpDir();
4040 TCollection_AsciiString aGenericName = (char*)aTmpDir.c_str();
4041 aGenericName += "NETGEN_";
4043 aGenericName += getpid();
4045 aGenericName += _getpid();
4047 aGenericName += "_";
4048 aGenericName += Abs((Standard_Integer)(long) aGenericName.ToCString());
4049 aGenericName += ".out";
4051 return aGenericName.ToCString();
4054 //================================================================================
4056 * \brief Remove "test.out" and "problemfaces" files in current directory
4058 //================================================================================
4060 void NETGENPlugin_NetgenLibWrapper::RemoveTmpFiles()
4062 bool rm = SMESH_File("test.out").remove() ;
4064 if ( rm && netgen::testout && instanceCounter() == 0 )
4066 delete netgen::testout;
4067 netgen::testout = 0;
4070 SMESH_File("problemfaces").remove();
4071 SMESH_File("occmesh.rep").remove();
4074 //================================================================================
4076 * \brief Remove file with netgen output
4078 //================================================================================
4080 void NETGENPlugin_NetgenLibWrapper::removeOutputFile()
4082 if ( !_outputFileName.empty() )
4086 delete netgen::mycout;
4087 netgen::mycout = _ngcout;
4088 netgen::myerr = _ngcerr;
4091 string tmpDir = SALOMEDS_Tool::GetDirFromPath ( _outputFileName );
4092 string aFileName = SALOMEDS_Tool::GetNameFromPath( _outputFileName ) + ".out";
4093 SALOMEDS::ListOfFileNames_var aFiles = new SALOMEDS::ListOfFileNames;
4095 aFiles[0] = aFileName.c_str();
4097 SALOMEDS_Tool::RemoveTemporaryFiles( tmpDir.c_str(), aFiles.in(), true );