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 #if defined(NETGEN_V5) && defined(WIN32)
88 extern MeshingParameters mparam;
89 #if defined(NETGEN_V5) && defined(WIN32)
92 extern volatile multithreadt multithread;
94 #if defined(NETGEN_V5) && defined(WIN32)
97 extern bool merge_solids;
99 // values used for occgeo.facemeshstatus
100 enum EFaceMeshStatus { FACE_NOT_TREATED = 0,
112 using namespace nglib;
116 #define nodeVec_ACCESS(index) ((SMDS_MeshNode*) nodeVec.at((index)))
118 #define nodeVec_ACCESS(index) ((SMDS_MeshNode*) nodeVec[index])
121 #define NGPOINT_COORDS(p) p(0),p(1),p(2)
124 // dump elements added to ng mesh
125 //#define DUMP_SEGMENTS
126 //#define DUMP_TRIANGLES
127 //#define DUMP_TRIANGLES_SCRIPT "/tmp/trias.py" //!< debug AddIntVerticesInSolids()
130 TopTools_IndexedMapOfShape ShapesWithLocalSize;
131 std::map<int,double> VertexId2LocalSize;
132 std::map<int,double> EdgeId2LocalSize;
133 std::map<int,double> FaceId2LocalSize;
134 std::map<int,double> SolidId2LocalSize;
136 std::vector<SMESHUtils::ControlPnt> ControlPoints;
137 std::set<int> ShapesWithControlPoints; // <-- allows calling SetLocalSize() several times w/o recomputing ControlPoints
139 //=============================================================================
143 //=============================================================================
145 NETGENPlugin_Mesher::NETGENPlugin_Mesher (SMESH_Mesh* mesh,
146 const TopoDS_Shape& aShape,
152 _fineness(NETGENPlugin_Hypothesis::GetDefaultFineness()),
153 _isViscousLayers2D(false),
162 SetDefaultParameters();
163 ShapesWithLocalSize.Clear();
164 VertexId2LocalSize.clear();
165 EdgeId2LocalSize.clear();
166 FaceId2LocalSize.clear();
167 SolidId2LocalSize.clear();
168 ControlPoints.clear();
169 ShapesWithControlPoints.clear();
172 //================================================================================
176 //================================================================================
178 NETGENPlugin_Mesher::~NETGENPlugin_Mesher()
186 //================================================================================
188 * Set pointer to NETGENPlugin_Mesher* field of the holder, that will be
189 * nullified at destruction of this
191 //================================================================================
193 void NETGENPlugin_Mesher::SetSelfPointer( NETGENPlugin_Mesher ** ptr )
204 //================================================================================
206 * \brief Initialize global NETGEN parameters with default values
208 //================================================================================
210 void NETGENPlugin_Mesher::SetDefaultParameters()
212 netgen::MeshingParameters& mparams = netgen::mparam;
213 // maximal mesh edge size
214 mparams.maxh = 0;//NETGENPlugin_Hypothesis::GetDefaultMaxSize();
216 // minimal number of segments per edge
217 mparams.segmentsperedge = NETGENPlugin_Hypothesis::GetDefaultNbSegPerEdge();
218 // rate of growth of size between elements
219 mparams.grading = NETGENPlugin_Hypothesis::GetDefaultGrowthRate();
220 // safety factor for curvatures (elements per radius)
221 mparams.curvaturesafety = NETGENPlugin_Hypothesis::GetDefaultNbSegPerRadius();
222 // create elements of second order
223 mparams.secondorder = NETGENPlugin_Hypothesis::GetDefaultSecondOrder();
224 // quad-dominated surface meshing
228 mparams.quad = NETGENPlugin_Hypothesis_2D::GetDefaultQuadAllowed();
229 _fineness = NETGENPlugin_Hypothesis::GetDefaultFineness();
230 mparams.uselocalh = NETGENPlugin_Hypothesis::GetDefaultSurfaceCurvature();
231 netgen::merge_solids = NETGENPlugin_Hypothesis::GetDefaultFuseEdges();
234 //=============================================================================
238 //=============================================================================
240 void SetLocalSize(TopoDS_Shape GeomShape, double LocalSize)
242 if ( GeomShape.IsNull() ) return;
243 TopAbs_ShapeEnum GeomType = GeomShape.ShapeType();
244 if (GeomType == TopAbs_COMPOUND) {
245 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()) {
246 SetLocalSize(it.Value(), LocalSize);
251 if (! ShapesWithLocalSize.Contains(GeomShape))
252 key = ShapesWithLocalSize.Add(GeomShape);
254 key = ShapesWithLocalSize.FindIndex(GeomShape);
255 if (GeomType == TopAbs_VERTEX) {
256 VertexId2LocalSize[key] = LocalSize;
257 } else if (GeomType == TopAbs_EDGE) {
258 EdgeId2LocalSize[key] = LocalSize;
259 } else if (GeomType == TopAbs_FACE) {
260 FaceId2LocalSize[key] = LocalSize;
261 } else if (GeomType == TopAbs_SOLID) {
262 SolidId2LocalSize[key] = LocalSize;
266 //=============================================================================
268 * Pass parameters to NETGEN
270 //=============================================================================
271 void NETGENPlugin_Mesher::SetParameters(const NETGENPlugin_Hypothesis* hyp)
275 netgen::MeshingParameters& mparams = netgen::mparam;
276 // Initialize global NETGEN parameters:
277 // maximal mesh segment size
278 mparams.maxh = hyp->GetMaxSize();
279 // maximal mesh element linear size
280 mparams.minh = hyp->GetMinSize();
281 // minimal number of segments per edge
282 mparams.segmentsperedge = hyp->GetNbSegPerEdge();
283 // rate of growth of size between elements
284 mparams.grading = hyp->GetGrowthRate();
285 // safety factor for curvatures (elements per radius)
286 mparams.curvaturesafety = hyp->GetNbSegPerRadius();
287 // create elements of second order
288 mparams.secondorder = hyp->GetSecondOrder() ? 1 : 0;
289 // quad-dominated surface meshing
290 mparams.quad = hyp->GetQuadAllowed() ? 1 : 0;
291 _optimize = hyp->GetOptimize();
292 _fineness = hyp->GetFineness();
293 mparams.uselocalh = hyp->GetSurfaceCurvature();
294 netgen::merge_solids = hyp->GetFuseEdges();
297 mparams.meshsizefilename= hyp->GetMeshSizeFile().empty() ? 0 : hyp->GetMeshSizeFile().c_str();
299 SMESH_Gen_i* smeshGen_i = SMESH_Gen_i::GetSMESHGen();
300 CORBA::Object_var anObject = smeshGen_i->GetNS()->Resolve("/myStudyManager");
301 SALOMEDS::StudyManager_var aStudyMgr = SALOMEDS::StudyManager::_narrow(anObject);
302 SALOMEDS::Study_var myStudy = aStudyMgr->GetStudyByID(hyp->GetStudyId());
303 if ( !myStudy->_is_nil() )
305 const NETGENPlugin_Hypothesis::TLocalSize localSizes = hyp->GetLocalSizesAndEntries();
306 NETGENPlugin_Hypothesis::TLocalSize::const_iterator it = localSizes.begin();
307 for ( ; it != localSizes.end() ; it++)
309 std::string entry = (*it).first;
310 double val = (*it).second;
312 GEOM::GEOM_Object_var aGeomObj;
313 SALOMEDS::SObject_var aSObj = myStudy->FindObjectID( entry.c_str() );
314 if ( !aSObj->_is_nil() ) {
315 CORBA::Object_var obj = aSObj->GetObject();
316 aGeomObj = GEOM::GEOM_Object::_narrow(obj);
319 TopoDS_Shape S = smeshGen_i->GeomObjectToShape( aGeomObj.in() );
320 ::SetLocalSize(S, val);
326 //=============================================================================
328 * Pass simple parameters to NETGEN
330 //=============================================================================
332 void NETGENPlugin_Mesher::SetParameters(const NETGENPlugin_SimpleHypothesis_2D* hyp)
336 SetDefaultParameters();
339 //=============================================================================
341 * Link - a pair of integer numbers
343 //=============================================================================
347 Link(int _n1, int _n2) : n1(_n1), n2(_n2) {}
348 Link() : n1(0), n2(0) {}
349 bool Contains( int n ) const { return n == n1 || n == n2; }
350 bool IsConnected( const Link& other ) const
352 return (( Contains( other.n1 ) || Contains( other.n2 )) && ( this != &other ));
356 int HashCode(const Link& aLink, int aLimit)
358 return HashCode(aLink.n1 + aLink.n2, aLimit);
361 Standard_Boolean IsEqual(const Link& aLink1, const Link& aLink2)
363 return (( aLink1.n1 == aLink2.n1 && aLink1.n2 == aLink2.n2 ) ||
364 ( aLink1.n1 == aLink2.n2 && aLink1.n2 == aLink2.n1 ));
369 //================================================================================
371 * \brief return id of netgen point corresponding to SMDS node
373 //================================================================================
374 typedef map< const SMDS_MeshNode*, int > TNode2IdMap;
376 int ngNodeId( const SMDS_MeshNode* node,
377 netgen::Mesh& ngMesh,
378 TNode2IdMap& nodeNgIdMap)
380 int newNgId = ngMesh.GetNP() + 1;
382 TNode2IdMap::iterator node_id = nodeNgIdMap.insert( make_pair( node, newNgId )).first;
384 if ( node_id->second == newNgId)
386 #if defined(DUMP_SEGMENTS) || defined(DUMP_TRIANGLES)
387 cout << "Ng " << newNgId << " - " << node;
389 netgen::MeshPoint p( netgen::Point<3> (node->X(), node->Y(), node->Z()) );
390 ngMesh.AddPoint( p );
392 return node_id->second;
395 //================================================================================
397 * \brief Return computed EDGEs connected to the given one
399 //================================================================================
401 list< TopoDS_Edge > getConnectedEdges( const TopoDS_Edge& edge,
402 const TopoDS_Face& face,
403 const set< SMESH_subMesh* > & computedSM,
404 const SMESH_MesherHelper& helper,
405 map< SMESH_subMesh*, set< int > >& addedEdgeSM2Faces)
408 list< TopoDS_Edge > edges;
409 list< int > nbEdgesInWire;
410 /*int nbWires =*/ SMESH_Block::GetOrderedEdges( face, edges, nbEdgesInWire);
412 // find <edge> within <edges>
413 list< TopoDS_Edge >::iterator eItFwd = edges.begin();
414 for ( ; eItFwd != edges.end(); ++eItFwd )
415 if ( edge.IsSame( *eItFwd ))
417 if ( eItFwd == edges.end()) return list< TopoDS_Edge>();
419 if ( eItFwd->Orientation() >= TopAbs_INTERNAL )
421 // connected INTERNAL edges returned from GetOrderedEdges() are wrongly oriented
422 // so treat each INTERNAL edge separately
423 TopoDS_Edge e = *eItFwd;
425 edges.push_back( e );
429 // get all computed EDGEs connected to <edge>
431 list< TopoDS_Edge >::iterator eItBack = eItFwd, ePrev;
432 TopoDS_Vertex vCommon;
433 TopTools_MapOfShape eAdded; // map used not to add a seam edge twice to <edges>
436 // put edges before <edge> to <edges> back
437 while ( edges.begin() != eItFwd )
438 edges.splice( edges.end(), edges, edges.begin() );
442 while ( ++eItFwd != edges.end() )
444 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( *eItFwd );
446 bool connected = TopExp::CommonVertex( *ePrev, *eItFwd, vCommon );
447 bool computed = sm->IsMeshComputed();
448 bool added = addedEdgeSM2Faces[ sm ].count( helper.GetSubShapeID() );
449 bool doubled = !eAdded.Add( *eItFwd );
450 bool orientOK = (( ePrev ->Orientation() < TopAbs_INTERNAL ) ==
451 ( eItFwd->Orientation() < TopAbs_INTERNAL ) );
452 if ( !connected || !computed || !orientOK || added || doubled )
454 // stop advancement; move edges from tail to head
455 while ( edges.back() != *ePrev )
456 edges.splice( edges.begin(), edges, --edges.end() );
462 while ( eItBack != edges.begin() )
466 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( *eItBack );
468 bool connected = TopExp::CommonVertex( *ePrev, *eItBack, vCommon );
469 bool computed = sm->IsMeshComputed();
470 bool added = addedEdgeSM2Faces[ sm ].count( helper.GetSubShapeID() );
471 bool doubled = !eAdded.Add( *eItBack );
472 bool orientOK = (( ePrev ->Orientation() < TopAbs_INTERNAL ) ==
473 ( eItBack->Orientation() < TopAbs_INTERNAL ) );
474 if ( !connected || !computed || !orientOK || added || doubled)
477 edges.erase( edges.begin(), ePrev );
481 if ( edges.front() != edges.back() )
483 // assure that the 1st vertex is meshed
484 TopoDS_Edge eLast = edges.back();
485 while ( !SMESH_Algo::VertexNode( SMESH_MesherHelper::IthVertex( 0, edges.front()), helper.GetMeshDS())
487 edges.front() != eLast )
488 edges.splice( edges.end(), edges, edges.begin() );
493 //================================================================================
495 * \brief Make triangulation of a shape precise enough
497 //================================================================================
499 void updateTriangulation( const TopoDS_Shape& shape )
501 // static set< Poly_Triangulation* > updated;
503 // TopLoc_Location loc;
504 // TopExp_Explorer fExp( shape, TopAbs_FACE );
505 // for ( ; fExp.More(); fExp.Next() )
507 // Handle(Poly_Triangulation) triangulation =
508 // BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
509 // if ( triangulation.IsNull() ||
510 // updated.insert( triangulation.operator->() ).second )
512 // BRepTools::Clean (shape);
515 BRepMesh_IncrementalMesh e(shape, 0.01, true);
517 catch (Standard_Failure)
520 // updated.erase( triangulation.operator->() );
521 // triangulation = BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
522 // updated.insert( triangulation.operator->() );
526 //================================================================================
528 * \brief Returns a medium node either existing in SMESH of created by NETGEN
529 * \param [in] corner1 - corner node 1
530 * \param [in] corner2 - corner node 2
531 * \param [in] defaultMedium - the node created by NETGEN
532 * \param [in] helper - holder of medium nodes existing in SMESH
533 * \return const SMDS_MeshNode* - the result node
535 //================================================================================
537 const SMDS_MeshNode* mediumNode( const SMDS_MeshNode* corner1,
538 const SMDS_MeshNode* corner2,
539 const SMDS_MeshNode* defaultMedium,
540 const SMESH_MesherHelper* helper)
544 TLinkNodeMap::const_iterator l2n =
545 helper->GetTLinkNodeMap().find( SMESH_TLink( corner1, corner2 ));
546 if ( l2n != helper->GetTLinkNodeMap().end() )
547 defaultMedium = l2n->second;
549 return defaultMedium;
552 //================================================================================
554 * \brief Assure that mesh on given shapes is quadratic
556 //================================================================================
558 // void makeQuadratic( const TopTools_IndexedMapOfShape& shapes,
559 // SMESH_Mesh* mesh )
561 // for ( int i = 1; i <= shapes.Extent(); ++i )
563 // SMESHDS_SubMesh* smDS = mesh->GetMeshDS()->MeshElements( shapes(i) );
564 // if ( !smDS ) continue;
565 // SMDS_ElemIteratorPtr elemIt = smDS->GetElements();
566 // if ( !elemIt->more() ) continue;
567 // const SMDS_MeshElement* e = elemIt->next();
568 // if ( !e || e->IsQuadratic() )
571 // TIDSortedElemSet elems;
572 // elems.insert( e );
573 // while ( elemIt->more() )
574 // elems.insert( elems.end(), elemIt->next() );
576 // SMESH_MeshEditor( mesh ).ConvertToQuadratic( /*3d=*/false, elems, /*biQuad=*/false );
580 //================================================================================
582 * \brief Restrict size of elements on the given edge
584 //================================================================================
586 void setLocalSize(const TopoDS_Edge& edge,
590 if ( size <= std::numeric_limits<double>::min() )
592 Standard_Real u1, u2;
593 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, u1, u2);
594 if ( curve.IsNull() )
596 TopoDS_Iterator vIt( edge );
597 if ( !vIt.More() ) return;
598 gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex( vIt.Value() ));
599 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
603 const int nb = (int)( 1.5 * SMESH_Algo::EdgeLength( edge ) / size );
604 Standard_Real delta = (u2-u1)/nb;
605 for(int i=0; i<nb; i++)
607 Standard_Real u = u1 + delta*i;
608 gp_Pnt p = curve->Value(u);
609 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
610 netgen::Point3d pi(p.X(), p.Y(), p.Z());
611 double resultSize = mesh.GetH(pi);
612 if ( resultSize - size > 0.1*size )
613 // netgen does restriction iff oldH/newH > 1.2 (localh.cpp:136)
614 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), resultSize/1.201 );
620 //================================================================================
622 * \brief Set local size on shapes defined by SetParameters()
624 //================================================================================
626 void NETGENPlugin_Mesher::SetLocalSize( netgen::OCCGeometry& occgeo,
627 netgen::Mesh& ngMesh )
629 for(std::map<int,double>::const_iterator it=EdgeId2LocalSize.begin(); it!=EdgeId2LocalSize.end(); it++)
631 int key = (*it).first;
632 double hi = (*it).second;
633 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
634 setLocalSize( TopoDS::Edge(shape), hi, ngMesh );
636 for(std::map<int,double>::const_iterator it=VertexId2LocalSize.begin(); it!=VertexId2LocalSize.end(); it++)
638 int key = (*it).first;
639 double hi = (*it).second;
640 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
641 gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex(shape) );
642 NETGENPlugin_Mesher::RestrictLocalSize( ngMesh, p.XYZ(), hi );
644 for(map<int,double>::const_iterator it=FaceId2LocalSize.begin(); it!=FaceId2LocalSize.end(); it++)
646 int key = (*it).first;
647 double val = (*it).second;
648 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
649 int faceNgID = occgeo.fmap.FindIndex(shape);
652 occgeo.SetFaceMaxH(faceNgID, val);
653 for ( TopExp_Explorer edgeExp( shape, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
654 setLocalSize( TopoDS::Edge( edgeExp.Current() ), val, ngMesh );
656 else if ( !ShapesWithControlPoints.count( key ))
658 SMESHUtils::createPointsSampleFromFace( TopoDS::Face( shape ), val, ControlPoints );
659 ShapesWithControlPoints.insert( key );
662 for(map<int,double>::const_iterator it=SolidId2LocalSize.begin(); it!=SolidId2LocalSize.end(); it++)
664 int key = (*it).first;
665 double val = (*it).second;
666 if ( !ShapesWithControlPoints.count( key ))
668 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
669 SMESHUtils::createPointsSampleFromSolid( TopoDS::Solid( shape ), val, ControlPoints );
670 ShapesWithControlPoints.insert( key );
674 if ( !ControlPoints.empty() )
676 for ( size_t i = 0; i < ControlPoints.size(); ++i )
677 NETGENPlugin_Mesher::RestrictLocalSize( ngMesh, ControlPoints[i].XYZ(), ControlPoints[i].Size() );
681 //================================================================================
683 * \brief Initialize netgen::OCCGeometry with OCCT shape
685 //================================================================================
687 void NETGENPlugin_Mesher::PrepareOCCgeometry(netgen::OCCGeometry& occgeo,
688 const TopoDS_Shape& shape,
690 list< SMESH_subMesh* > * meshedSM,
691 NETGENPlugin_Internals* intern)
693 updateTriangulation( shape );
696 BRepBndLib::Add (shape, bb);
697 double x1,y1,z1,x2,y2,z2;
698 bb.Get (x1,y1,z1,x2,y2,z2);
699 netgen::Point<3> p1 = netgen::Point<3> (x1,y1,z1);
700 netgen::Point<3> p2 = netgen::Point<3> (x2,y2,z2);
701 occgeo.boundingbox = netgen::Box<3> (p1,p2);
703 occgeo.shape = shape;
706 // fill maps of shapes of occgeo with not yet meshed subshapes
708 // get root submeshes
709 list< SMESH_subMesh* > rootSM;
710 const int shapeID = mesh.GetMeshDS()->ShapeToIndex( shape );
711 if ( shapeID > 0 ) { // SMESH_subMesh with ID 0 may exist, don't use it!
712 rootSM.push_back( mesh.GetSubMesh( shape ));
715 for ( TopoDS_Iterator it( shape ); it.More(); it.Next() )
716 rootSM.push_back( mesh.GetSubMesh( it.Value() ));
721 // add subshapes of empty submeshes
722 list< SMESH_subMesh* >::iterator rootIt = rootSM.begin(), rootEnd = rootSM.end();
723 for ( ; rootIt != rootEnd; ++rootIt ) {
724 SMESH_subMesh * root = *rootIt;
725 SMESH_subMeshIteratorPtr smIt = root->getDependsOnIterator(/*includeSelf=*/true,
726 /*complexShapeFirst=*/true);
727 // to find a right orientation of subshapes (PAL20462)
728 TopTools_IndexedMapOfShape subShapes;
729 TopExp::MapShapes(root->GetSubShape(), subShapes);
730 while ( smIt->more() )
732 SMESH_subMesh* sm = smIt->next();
733 TopoDS_Shape shape = sm->GetSubShape();
734 totNbFaces += ( shape.ShapeType() == TopAbs_FACE );
735 if ( intern && intern->isShapeToPrecompute( shape ))
737 if ( !meshedSM || sm->IsEmpty() )
739 if ( shape.ShapeType() != TopAbs_VERTEX )
740 shape = subShapes( subShapes.FindIndex( shape ));// shape -> index -> oriented shape
741 if ( shape.Orientation() >= TopAbs_INTERNAL )
742 shape.Orientation( TopAbs_FORWARD ); // isuue 0020676
743 switch ( shape.ShapeType() ) {
744 case TopAbs_FACE : occgeo.fmap.Add( shape ); break;
745 case TopAbs_EDGE : occgeo.emap.Add( shape ); break;
746 case TopAbs_VERTEX: occgeo.vmap.Add( shape ); break;
747 case TopAbs_SOLID :occgeo.somap.Add( shape ); break;
751 // collect submeshes of meshed shapes
754 const int dim = SMESH_Gen::GetShapeDim( shape );
755 meshedSM[ dim ].push_back( sm );
759 occgeo.facemeshstatus.SetSize (totNbFaces);
760 occgeo.facemeshstatus = 0;
761 occgeo.face_maxh_modified.SetSize(totNbFaces);
762 occgeo.face_maxh_modified = 0;
763 occgeo.face_maxh.SetSize(totNbFaces);
764 occgeo.face_maxh = netgen::mparam.maxh;
767 //================================================================================
769 * \brief Return a default min size value suitable for the given geometry.
771 //================================================================================
773 double NETGENPlugin_Mesher::GetDefaultMinSize(const TopoDS_Shape& geom,
774 const double maxSize)
776 updateTriangulation( geom );
780 const int* pi[4] = { &i1, &i2, &i3, &i1 };
783 TopExp_Explorer fExp( geom, TopAbs_FACE );
784 for ( ; fExp.More(); fExp.Next() )
786 Handle(Poly_Triangulation) triangulation =
787 BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
788 if ( triangulation.IsNull() ) continue;
789 const double fTol = BRep_Tool::Tolerance( TopoDS::Face( fExp.Current() ));
790 const TColgp_Array1OfPnt& points = triangulation->Nodes();
791 const Poly_Array1OfTriangle& trias = triangulation->Triangles();
792 for ( int iT = trias.Lower(); iT <= trias.Upper(); ++iT )
794 trias(iT).Get( i1, i2, i3 );
795 for ( int j = 0; j < 3; ++j )
797 double dist2 = points(*pi[j]).SquareDistance( points( *pi[j+1] ));
798 if ( dist2 < minh && fTol*fTol < dist2 )
800 bb.Add( points(*pi[j]));
804 if ( minh > 0.25 * bb.SquareExtent() ) // simple geometry, rough triangulation
806 minh = 1e-3 * sqrt( bb.SquareExtent());
807 //cout << "BND BOX minh = " <<minh << endl;
811 minh = sqrt( minh ); // triangulation for visualization is rather fine
812 //cout << "TRIANGULATION minh = " <<minh << endl;
814 if ( minh > 0.5 * maxSize )
820 //================================================================================
822 * \brief Restrict size of elements at a given point
824 //================================================================================
826 void NETGENPlugin_Mesher::RestrictLocalSize(netgen::Mesh& ngMesh,
829 const bool overrideMinH)
831 if ( size <= std::numeric_limits<double>::min() )
833 if ( netgen::mparam.minh > size )
837 ngMesh.SetMinimalH( size );
838 netgen::mparam.minh = size;
842 size = netgen::mparam.minh;
845 netgen::Point3d pi(p.X(), p.Y(), p.Z());
846 ngMesh.RestrictLocalH( pi, size );
849 //================================================================================
851 * \brief fill ngMesh with nodes and elements of computed submeshes
853 //================================================================================
855 bool NETGENPlugin_Mesher::FillNgMesh(netgen::OCCGeometry& occgeom,
856 netgen::Mesh& ngMesh,
857 vector<const SMDS_MeshNode*>& nodeVec,
858 const list< SMESH_subMesh* > & meshedSM,
859 SMESH_MesherHelper* quadHelper,
860 SMESH_ProxyMesh::Ptr proxyMesh)
862 TNode2IdMap nodeNgIdMap;
863 for ( size_t i = 1; i < nodeVec.size(); ++i )
864 nodeNgIdMap.insert( make_pair( nodeVec[i], i ));
866 TopTools_MapOfShape visitedShapes;
867 map< SMESH_subMesh*, set< int > > visitedEdgeSM2Faces;
868 set< SMESH_subMesh* > computedSM( meshedSM.begin(), meshedSM.end() );
870 SMESH_MesherHelper helper (*_mesh);
872 int faceNgID = ngMesh.GetNFD();
874 list< SMESH_subMesh* >::const_iterator smIt, smEnd = meshedSM.end();
875 for ( smIt = meshedSM.begin(); smIt != smEnd; ++smIt )
877 SMESH_subMesh* sm = *smIt;
878 if ( !visitedShapes.Add( sm->GetSubShape() ))
881 const SMESHDS_SubMesh * smDS = sm->GetSubMeshDS();
882 if ( !smDS ) continue;
884 switch ( sm->GetSubShape().ShapeType() )
886 case TopAbs_EDGE: { // EDGE
887 // ----------------------
888 TopoDS_Edge geomEdge = TopoDS::Edge( sm->GetSubShape() );
889 if ( geomEdge.Orientation() >= TopAbs_INTERNAL )
890 geomEdge.Orientation( TopAbs_FORWARD ); // issue 0020676
892 // Add ng segments for each not meshed FACE the EDGE bounds
893 PShapeIteratorPtr fIt = helper.GetAncestors( geomEdge, *sm->GetFather(), TopAbs_FACE );
894 while ( const TopoDS_Shape * anc = fIt->next() )
896 faceNgID = occgeom.fmap.FindIndex( *anc );
898 continue; // meshed face
900 int faceSMDSId = helper.GetMeshDS()->ShapeToIndex( *anc );
901 if ( visitedEdgeSM2Faces[ sm ].count( faceSMDSId ))
902 continue; // already treated EDGE
904 TopoDS_Face face = TopoDS::Face( occgeom.fmap( faceNgID ));
905 if ( face.Orientation() >= TopAbs_INTERNAL )
906 face.Orientation( TopAbs_FORWARD ); // issue 0020676
908 // get all meshed EDGEs of the FACE connected to geomEdge (issue 0021140)
909 helper.SetSubShape( face );
910 list< TopoDS_Edge > edges = getConnectedEdges( geomEdge, face, computedSM, helper,
911 visitedEdgeSM2Faces );
913 continue; // wrong ancestor?
915 // find out orientation of <edges> within <face>
916 TopoDS_Edge eNotSeam = edges.front();
917 if ( helper.HasSeam() )
919 list< TopoDS_Edge >::iterator eIt = edges.begin();
920 while ( helper.IsRealSeam( *eIt )) ++eIt;
921 if ( eIt != edges.end() )
924 TopAbs_Orientation fOri = helper.GetSubShapeOri( face, eNotSeam );
925 bool isForwad = ( fOri == eNotSeam.Orientation() || fOri >= TopAbs_INTERNAL );
927 // get all nodes from connected <edges>
928 const bool isQuad = smDS->IsQuadratic();
929 StdMeshers_FaceSide fSide( face, edges, _mesh, isForwad, isQuad );
930 const vector<UVPtStruct>& points = fSide.GetUVPtStruct();
931 if ( points.empty() )
932 return false; // invalid node params?
933 int i, nbSeg = fSide.NbSegments();
935 // remember EDGEs of fSide to treat only once
936 for ( int iE = 0; iE < fSide.NbEdges(); ++iE )
937 visitedEdgeSM2Faces[ helper.GetMesh()->GetSubMesh( fSide.Edge(iE )) ].insert(faceSMDSId);
939 double otherSeamParam = 0;
944 int prevNgId = ngNodeId( points[0].node, ngMesh, nodeNgIdMap );
946 for ( i = 0; i < nbSeg; ++i )
948 const UVPtStruct& p1 = points[ i ];
949 const UVPtStruct& p2 = points[ i+1 ];
951 if ( p1.node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX ) //an EDGE begins
954 if ( helper.IsRealSeam( p1.node->getshapeId() ))
956 TopoDS_Edge e = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
957 isSeam = helper.IsRealSeam( e );
960 otherSeamParam = helper.GetOtherParam( helper.GetPeriodicIndex() & 1 ? p2.u : p2.v );
967 seg[1] = prevNgId = ngNodeId( p2.node, ngMesh, nodeNgIdMap );
968 // node param on curve
969 seg.epgeominfo[ 0 ].dist = p1.param;
970 seg.epgeominfo[ 1 ].dist = p2.param;
972 seg.epgeominfo[ 0 ].u = p1.u;
973 seg.epgeominfo[ 0 ].v = p1.v;
974 seg.epgeominfo[ 1 ].u = p2.u;
975 seg.epgeominfo[ 1 ].v = p2.v;
977 //geomEdge = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
978 //seg.epgeominfo[ 0 ].edgenr = seg.epgeominfo[ 1 ].edgenr = occgeom.emap.FindIndex( geomEdge );
980 //seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
981 seg.si = faceNgID; // = geom.fmap.FindIndex (face);
982 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
983 ngMesh.AddSegment (seg);
985 SMESH_TNodeXYZ np1( p1.node ), np2( p2.node );
986 RestrictLocalSize( ngMesh, 0.5*(np1+np2), (np1-np2).Modulus() );
989 cout << "Segment: " << seg.edgenr << " on SMESH face " << helper.GetMeshDS()->ShapeToIndex( face ) << endl
990 << "\tface index: " << seg.si << endl
991 << "\tp1: " << seg[0] << endl
992 << "\tp2: " << seg[1] << endl
993 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
994 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
995 //<< "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
996 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
997 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl;
998 //<< "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
1002 if ( helper.GetPeriodicIndex() && 1 ) {
1003 seg.epgeominfo[ 0 ].u = otherSeamParam;
1004 seg.epgeominfo[ 1 ].u = otherSeamParam;
1005 swap (seg.epgeominfo[0].v, seg.epgeominfo[1].v);
1007 seg.epgeominfo[ 0 ].v = otherSeamParam;
1008 seg.epgeominfo[ 1 ].v = otherSeamParam;
1009 swap (seg.epgeominfo[0].u, seg.epgeominfo[1].u);
1011 swap( seg[0], seg[1] );
1012 swap( seg.epgeominfo[0].dist, seg.epgeominfo[1].dist );
1013 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1014 ngMesh.AddSegment( seg );
1015 #ifdef DUMP_SEGMENTS
1016 cout << "Segment: " << seg.edgenr << endl
1017 << "\t is SEAM (reverse) of the previous. "
1018 << " Other " << (helper.GetPeriodicIndex() && 1 ? "U" : "V")
1019 << " = " << otherSeamParam << endl;
1022 else if ( fOri == TopAbs_INTERNAL )
1024 swap( seg[0], seg[1] );
1025 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1026 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1027 ngMesh.AddSegment( seg );
1028 #ifdef DUMP_SEGMENTS
1029 cout << "Segment: " << seg.edgenr << endl << "\t is REVERSE of the previous" << endl;
1033 } // loop on geomEdge ancestors
1035 if ( quadHelper ) // remember medium nodes of sub-meshes
1037 SMDS_ElemIteratorPtr edges = smDS->GetElements();
1038 while ( edges->more() )
1040 const SMDS_MeshElement* e = edges->next();
1041 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshEdge*>( e )))
1047 } // case TopAbs_EDGE
1049 case TopAbs_FACE: { // FACE
1050 // ----------------------
1051 const TopoDS_Face& geomFace = TopoDS::Face( sm->GetSubShape() );
1052 helper.SetSubShape( geomFace );
1053 bool isInternalFace = ( geomFace.Orientation() == TopAbs_INTERNAL );
1055 // Find solids the geomFace bounds
1056 int solidID1 = 0, solidID2 = 0;
1057 StdMeshers_QuadToTriaAdaptor* quadAdaptor =
1058 dynamic_cast<StdMeshers_QuadToTriaAdaptor*>( proxyMesh.get() );
1061 solidID1 = occgeom.somap.FindIndex( quadAdaptor->GetShape() );
1065 PShapeIteratorPtr solidIt = helper.GetAncestors( geomFace, *sm->GetFather(), TopAbs_SOLID);
1066 while ( const TopoDS_Shape * solid = solidIt->next() )
1068 int id = occgeom.somap.FindIndex ( *solid );
1069 if ( solidID1 && id != solidID1 ) solidID2 = id;
1073 // Add ng face descriptors of meshed faces
1075 ngMesh.AddFaceDescriptor( netgen::FaceDescriptor( faceNgID, solidID1, solidID2, 0 ));
1077 // if second oreder is required, even already meshed faces must be passed to NETGEN
1078 int fID = occgeom.fmap.Add( geomFace );
1079 if ( occgeom.facemeshstatus.Size() < fID ) occgeom.facemeshstatus.SetSize( fID );
1080 occgeom.facemeshstatus[ fID-1 ] = netgen::FACE_MESHED_OK;
1081 while ( fID < faceNgID ) // geomFace is already in occgeom.fmap, add a copy
1083 fID = occgeom.fmap.Add( BRepBuilderAPI_Copy( geomFace, /*copyGeom=*/false ));
1084 if ( occgeom.facemeshstatus.Size() < fID ) occgeom.facemeshstatus.SetSize( fID );
1085 occgeom.facemeshstatus[ fID-1 ] = netgen::FACE_MESHED_OK;
1087 // Problem with the second order in a quadrangular mesh remains.
1088 // 1) All quadrangles generated by NETGEN are moved to an inexistent face
1089 // by FillSMesh() (find "AddFaceDescriptor")
1090 // 2) Temporary triangles generated by StdMeshers_QuadToTriaAdaptor
1091 // are on faces where quadrangles were.
1092 // Due to these 2 points, wrong geom faces are used while conversion to quadratic
1093 // of the mentioned above quadrangles and triangles
1095 // Orient the face correctly in solidID1 (issue 0020206)
1096 bool reverse = false;
1098 TopoDS_Shape solid = occgeom.somap( solidID1 );
1099 TopAbs_Orientation faceOriInSolid = helper.GetSubShapeOri( solid, geomFace );
1100 if ( faceOriInSolid >= 0 )
1102 helper.IsReversedSubMesh( TopoDS::Face( geomFace.Oriented( faceOriInSolid )));
1105 // Add surface elements
1107 netgen::Element2d tri(3);
1108 tri.SetIndex( faceNgID );
1109 SMESH_TNodeXYZ xyz[3];
1111 #ifdef DUMP_TRIANGLES
1112 cout << "SMESH face " << helper.GetMeshDS()->ShapeToIndex( geomFace )
1113 << " internal="<<isInternalFace << endl;
1116 smDS = proxyMesh->GetSubMesh( geomFace );
1118 SMDS_ElemIteratorPtr faces = smDS->GetElements();
1119 while ( faces->more() )
1121 const SMDS_MeshElement* f = faces->next();
1122 if ( f->NbNodes() % 3 != 0 ) // not triangle
1124 PShapeIteratorPtr solidIt=helper.GetAncestors(geomFace,*sm->GetFather(),TopAbs_SOLID);
1125 if ( const TopoDS_Shape * solid = solidIt->next() )
1126 sm = _mesh->GetSubMesh( *solid );
1127 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1128 smError.reset( new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH,"Not triangle sub-mesh"));
1129 smError->myBadElements.push_back( f );
1133 for ( int i = 0; i < 3; ++i )
1135 const SMDS_MeshNode* node = f->GetNode( i ), * inFaceNode=0;
1138 // get node UV on face
1139 int shapeID = node->getshapeId();
1140 if ( helper.IsSeamShape( shapeID ))
1142 if ( helper.IsSeamShape( f->GetNodeWrap( i+1 )->getshapeId() ))
1143 inFaceNode = f->GetNodeWrap( i-1 );
1145 inFaceNode = f->GetNodeWrap( i+1 );
1147 gp_XY uv = helper.GetNodeUV( geomFace, node, inFaceNode );
1149 int ind = reverse ? 3-i : i+1;
1150 tri.GeomInfoPi(ind).u = uv.X();
1151 tri.GeomInfoPi(ind).v = uv.Y();
1152 tri.PNum (ind) = ngNodeId( node, ngMesh, nodeNgIdMap );
1155 // pass a triangle size to NG size-map
1156 double size = ( ( xyz[0] - xyz[1] ).Modulus() +
1157 ( xyz[1] - xyz[2] ).Modulus() +
1158 ( xyz[2] - xyz[0] ).Modulus() ) / 3;
1159 gp_XYZ gc = ( xyz[0] + xyz[1] + xyz[2] ) / 3;
1160 RestrictLocalSize( ngMesh, gc, size, /*overrideMinH=*/false );
1162 ngMesh.AddSurfaceElement (tri);
1163 #ifdef DUMP_TRIANGLES
1164 cout << tri << endl;
1167 if ( isInternalFace )
1169 swap( tri[1], tri[2] );
1170 ngMesh.AddSurfaceElement (tri);
1171 #ifdef DUMP_TRIANGLES
1172 cout << tri << endl;
1177 if ( quadHelper ) // remember medium nodes of sub-meshes
1179 SMDS_ElemIteratorPtr faces = smDS->GetElements();
1180 while ( faces->more() )
1182 const SMDS_MeshElement* f = faces->next();
1183 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshFace*>( f )))
1189 } // case TopAbs_FACE
1191 case TopAbs_VERTEX: { // VERTEX
1192 // --------------------------
1193 // issue 0021405. Add node only if a VERTEX is shared by a not meshed EDGE,
1194 // else netgen removes a free node and nodeVector becomes invalid
1195 PShapeIteratorPtr ansIt = helper.GetAncestors( sm->GetSubShape(),
1199 while ( const TopoDS_Shape* e = ansIt->next() )
1201 SMESH_subMesh* eSub = helper.GetMesh()->GetSubMesh( *e );
1202 if (( toAdd = ( eSub->IsEmpty() && !SMESH_Algo::isDegenerated( TopoDS::Edge( *e )))))
1207 SMDS_NodeIteratorPtr nodeIt = smDS->GetNodes();
1208 if ( nodeIt->more() )
1209 ngNodeId( nodeIt->next(), ngMesh, nodeNgIdMap );
1215 } // loop on submeshes
1218 nodeVec.resize( ngMesh.GetNP() + 1 );
1219 TNode2IdMap::iterator node_NgId, nodeNgIdEnd = nodeNgIdMap.end();
1220 for ( node_NgId = nodeNgIdMap.begin(); node_NgId != nodeNgIdEnd; ++node_NgId)
1221 nodeVec[ node_NgId->second ] = node_NgId->first;
1226 //================================================================================
1228 * \brief Duplicate mesh faces on internal geom faces
1230 //================================================================================
1232 void NETGENPlugin_Mesher::FixIntFaces(const netgen::OCCGeometry& occgeom,
1233 netgen::Mesh& ngMesh,
1234 NETGENPlugin_Internals& internalShapes)
1236 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1238 // find ng indices of internal faces
1240 for ( int ngFaceID = 1; ngFaceID <= occgeom.fmap.Extent(); ++ngFaceID )
1242 int smeshID = meshDS->ShapeToIndex( occgeom.fmap( ngFaceID ));
1243 if ( internalShapes.isInternalShape( smeshID ))
1244 ngFaceIds.insert( ngFaceID );
1246 if ( !ngFaceIds.empty() )
1249 int i, nbFaces = ngMesh.GetNSE();
1250 for ( i = 1; i <= nbFaces; ++i)
1252 netgen::Element2d elem = ngMesh.SurfaceElement(i);
1253 if ( ngFaceIds.count( elem.GetIndex() ))
1255 swap( elem[1], elem[2] );
1256 ngMesh.AddSurfaceElement (elem);
1262 //================================================================================
1264 * \brief Tries to heal the mesh on a FACE. The FACE is supposed to be partially
1265 * meshed due to NETGEN failure
1266 * \param [in] occgeom - geometry
1267 * \param [in,out] ngMesh - the mesh to fix
1268 * \param [inout] faceID - ID of the FACE to fix the mesh on
1269 * \return bool - is mesh is or becomes OK
1271 //================================================================================
1273 bool NETGENPlugin_Mesher::FixFaceMesh(const netgen::OCCGeometry& occgeom,
1274 netgen::Mesh& ngMesh,
1277 // we address a case where the FACE is almost fully meshed except small holes
1278 // of usually triangular shape at FACE boundary (IPAL52861)
1280 // The case appeared to be not simple: holes only look triangular but
1281 // indeed are a self intersecting polygon. A reason of the bug was in coincident
1282 // NG points on a seam edge. But the code below is very nice, leave it for
1287 if ( occgeom.fmap.Extent() < faceID )
1289 //const TopoDS_Face& face = TopoDS::Face( occgeom.fmap( faceID ));
1291 // find free links on the FACE
1292 NCollection_Map<Link> linkMap;
1293 for ( int iF = 1; iF <= ngMesh.GetNSE(); ++iF )
1295 const netgen::Element2d& elem = ngMesh.SurfaceElement(iF);
1296 if ( faceID != elem.GetIndex() )
1298 int n0 = elem[ elem.GetNP() - 1 ];
1299 for ( int i = 0; i < elem.GetNP(); ++i )
1302 Link link( n0, n1 );
1303 if ( !linkMap.Add( link ))
1304 linkMap.Remove( link );
1308 // add/remove boundary links
1309 for ( int iSeg = 1; iSeg <= ngMesh.GetNSeg(); ++iSeg )
1311 const netgen::Segment& seg = ngMesh.LineSegment( iSeg );
1312 if ( seg.si != faceID ) // !edgeIDs.Contains( seg.edgenr ))
1314 Link link( seg[1], seg[0] ); // reverse!!!
1315 if ( !linkMap.Add( link ))
1316 linkMap.Remove( link );
1318 if ( linkMap.IsEmpty() )
1320 if ( linkMap.Extent() < 3 )
1323 // make triangles of the links
1325 netgen::Element2d tri(3);
1326 tri.SetIndex ( faceID );
1328 NCollection_Map<Link>::Iterator linkIt( linkMap );
1329 Link link1 = linkIt.Value();
1330 // look for a link connected to link1
1331 NCollection_Map<Link>::Iterator linkIt2 = linkIt;
1332 for ( linkIt2.Next(); linkIt2.More(); linkIt2.Next() )
1334 const Link& link2 = linkIt2.Value();
1335 if ( link2.IsConnected( link1 ))
1337 // look for a link connected to both link1 and link2
1338 NCollection_Map<Link>::Iterator linkIt3 = linkIt2;
1339 for ( linkIt3.Next(); linkIt3.More(); linkIt3.Next() )
1341 const Link& link3 = linkIt3.Value();
1342 if ( link3.IsConnected( link1 ) &&
1343 link3.IsConnected( link2 ) )
1348 tri[2] = ( link2.Contains( link1.n1 ) ? link2.n1 : link3.n1 );
1349 if ( tri[0] == tri[2] || tri[1] == tri[2] )
1351 ngMesh.AddSurfaceElement( tri );
1353 // prepare for the next tria search
1354 if ( linkMap.Extent() == 3 )
1356 linkMap.Remove( link3 );
1357 linkMap.Remove( link2 );
1359 linkMap.Remove( link1 );
1360 link1 = linkIt.Value();
1373 //================================================================================
1374 // define gp_XY_Subtracted pointer to function calling gp_XY::Subtracted(gp_XY)
1375 gp_XY_FunPtr(Subtracted);
1376 //gp_XY_FunPtr(Added);
1378 //================================================================================
1380 * \brief Evaluate distance between two 2d points along the surface
1382 //================================================================================
1384 double evalDist( const gp_XY& uv1,
1386 const Handle(Geom_Surface)& surf,
1387 const int stopHandler=-1)
1389 if ( stopHandler > 0 ) // continue recursion
1391 gp_XY mid = SMESH_MesherHelper::GetMiddleUV( surf, uv1, uv2 );
1392 return evalDist( uv1,mid, surf, stopHandler-1 ) + evalDist( mid,uv2, surf, stopHandler-1 );
1394 double dist3D = surf->Value( uv1.X(), uv1.Y() ).Distance( surf->Value( uv2.X(), uv2.Y() ));
1395 if ( stopHandler == 0 ) // stop recursion
1398 // start recursion if necessary
1399 double dist2D = SMESH_MesherHelper::ApplyIn2D(surf, uv1, uv2, gp_XY_Subtracted, 0).Modulus();
1400 if ( fabs( dist3D - dist2D ) < dist2D * 1e-10 )
1401 return dist3D; // equal parametrization of a planar surface
1403 return evalDist( uv1, uv2, surf, 3 ); // start recursion
1406 //================================================================================
1408 * \brief Data of vertex internal in geom face
1410 //================================================================================
1414 gp_XY uv; //!< UV in face parametric space
1415 int ngId; //!< ng id of corrsponding node
1416 gp_XY uvClose; //!< UV of closest boundary node
1417 int ngIdClose; //!< ng id of closest boundary node
1420 //================================================================================
1422 * \brief Data of vertex internal in solid
1424 //================================================================================
1428 int ngId; //!< ng id of corresponding node
1429 int ngIdClose; //!< ng id of closest 2d mesh element
1430 int ngIdCloseN; //!< ng id of closest node of the closest 2d mesh element
1433 inline double dist2(const netgen::MeshPoint& p1, const netgen::MeshPoint& p2)
1435 return gp_Pnt( NGPOINT_COORDS(p1)).SquareDistance( gp_Pnt( NGPOINT_COORDS(p2)));
1439 //================================================================================
1441 * \brief Make netgen take internal vertices in faces into account by adding
1442 * segments including internal vertices
1444 * This function works in supposition that 1D mesh is already computed in ngMesh
1446 //================================================================================
1448 void NETGENPlugin_Mesher::AddIntVerticesInFaces(const netgen::OCCGeometry& occgeom,
1449 netgen::Mesh& ngMesh,
1450 vector<const SMDS_MeshNode*>& nodeVec,
1451 NETGENPlugin_Internals& internalShapes)
1453 if ((int) nodeVec.size() < ngMesh.GetNP() )
1454 nodeVec.resize( ngMesh.GetNP(), 0 );
1456 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1457 SMESH_MesherHelper helper( internalShapes.getMesh() );
1459 const map<int,list<int> >& face2Vert = internalShapes.getFacesWithVertices();
1460 map<int,list<int> >::const_iterator f2v = face2Vert.begin();
1461 for ( ; f2v != face2Vert.end(); ++f2v )
1463 const TopoDS_Face& face = TopoDS::Face( meshDS->IndexToShape( f2v->first ));
1464 if ( face.IsNull() ) continue;
1465 int faceNgID = occgeom.fmap.FindIndex (face);
1466 if ( faceNgID < 0 ) continue;
1468 TopLoc_Location loc;
1469 Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc);
1471 helper.SetSubShape( face );
1472 helper.SetElementsOnShape( true );
1474 // Get data of internal vertices and add them to ngMesh
1476 multimap< double, TIntVData > dist2VData; // sort vertices by distance from boundary nodes
1478 int i, nbSegInit = ngMesh.GetNSeg();
1480 // boundary characteristics
1481 double totSegLen2D = 0;
1484 const list<int>& iVertices = f2v->second;
1485 list<int>::const_iterator iv = iVertices.begin();
1486 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1489 // get node on vertex
1490 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1491 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1494 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1495 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1496 nV = SMESH_Algo::VertexNode( V, meshDS );
1497 if ( !nV ) continue;
1500 netgen::MeshPoint mp( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1501 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1502 vData.ngId = ngMesh.GetNP();
1503 nodeVec.push_back( nV );
1507 vData.uv = helper.GetNodeUV( face, nV, 0, &uvOK );
1508 if ( !uvOK ) helper.CheckNodeUV( face, nV, vData.uv, BRep_Tool::Tolerance(V),/*force=*/1);
1510 // loop on all segments of the face to find the node closest to vertex and to count
1511 // average segment 2d length
1512 double closeDist2 = numeric_limits<double>::max(), dist2;
1514 for (i = 1; i <= ngMesh.GetNSeg(); ++i)
1516 netgen::Segment & seg = ngMesh.LineSegment(i);
1517 if ( seg.si != faceNgID ) continue;
1519 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1521 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1522 if ( ngIdLast == seg[ iEnd ] ) continue;
1523 dist2 = helper.ApplyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1524 if ( dist2 < closeDist2 )
1525 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1526 ngIdLast = seg[ iEnd ];
1530 totSegLen2D += helper.ApplyIn2D(surf, uv[0], uv[1], gp_XY_Subtracted, false).Modulus();
1534 dist2VData.insert( make_pair( closeDist2, vData ));
1537 if ( totNbSeg == 0 ) break;
1538 double avgSegLen2d = totSegLen2D / totNbSeg;
1540 // Loop on vertices to add segments
1542 multimap< double, TIntVData >::iterator dist_vData = dist2VData.begin();
1543 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1545 double closeDist2 = dist_vData->first, dist2;
1546 TIntVData & vData = dist_vData->second;
1548 // try to find more close node among segments added for internal vertices
1549 for (i = nbSegInit+1; i <= ngMesh.GetNSeg(); ++i)
1551 netgen::Segment & seg = ngMesh.LineSegment(i);
1552 if ( seg.si != faceNgID ) continue;
1554 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1556 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1557 dist2 = helper.ApplyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1558 if ( dist2 < closeDist2 )
1559 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1562 // decide whether to use the closest node as the second end of segment or to
1563 // create a new point
1564 int segEnd1 = vData.ngId;
1565 int segEnd2 = vData.ngIdClose; // to use closest node
1566 gp_XY uvV = vData.uv, uvP = vData.uvClose;
1567 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1568 double nodeDist2D = sqrt( closeDist2 );
1569 double nodeDist3D = evalDist( vData.uv, vData.uvClose, surf );
1570 bool avgLenOK = ( avgSegLen2d < 0.75 * nodeDist2D );
1571 bool hintLenOK = ( segLenHint < 0.75 * nodeDist3D );
1572 //cout << "uvV " << uvV.X() <<","<<uvV.Y() << " ";
1573 if ( hintLenOK || avgLenOK )
1575 // create a point between the closest node and V
1578 double r = min( 0.5, ( hintLenOK ? segLenHint/nodeDist3D : avgSegLen2d/nodeDist2D ));
1579 // direction from V to closet node in 2D
1580 gp_Dir2d v2n( helper.ApplyIn2D(surf, uvP, uvV, gp_XY_Subtracted, false ));
1582 uvP = vData.uv + r * nodeDist2D * v2n.XY();
1583 gp_Pnt P = surf->Value( uvP.X(), uvP.Y() ).Transformed( loc );
1585 netgen::MeshPoint mp( netgen::Point<3> (P.X(), P.Y(), P.Z()));
1586 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1587 segEnd2 = ngMesh.GetNP();
1588 //cout << "Middle " << r << " uv " << uvP.X() << "," << uvP.Y() << "( " << ngMesh.Point(segEnd2).X()<<","<<ngMesh.Point(segEnd2).Y()<<","<<ngMesh.Point(segEnd2).Z()<<" )"<< endl;
1589 SMDS_MeshNode * nP = helper.AddNode(P.X(), P.Y(), P.Z());
1590 nodeVec.push_back( nP );
1592 //else cout << "at Node " << " uv " << uvP.X() << "," << uvP.Y() << endl;
1595 netgen::Segment seg;
1597 if ( segEnd1 > segEnd2 ) swap( segEnd1, segEnd2 ), swap( uvV, uvP );
1598 seg[0] = segEnd1; // ng node id
1599 seg[1] = segEnd2; // ng node id
1600 seg.edgenr = ngMesh.GetNSeg() + 1;// segment id
1603 seg.epgeominfo[ 0 ].dist = 0; // param on curve
1604 seg.epgeominfo[ 0 ].u = uvV.X();
1605 seg.epgeominfo[ 0 ].v = uvV.Y();
1606 seg.epgeominfo[ 1 ].dist = 1; // param on curve
1607 seg.epgeominfo[ 1 ].u = uvP.X();
1608 seg.epgeominfo[ 1 ].v = uvP.Y();
1610 // seg.epgeominfo[ 0 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1611 // seg.epgeominfo[ 1 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1613 ngMesh.AddSegment (seg);
1615 // add reverse segment
1616 swap( seg[0], seg[1] );
1617 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1618 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1619 ngMesh.AddSegment (seg);
1625 //================================================================================
1627 * \brief Make netgen take internal vertices in solids into account by adding
1628 * faces including internal vertices
1630 * This function works in supposition that 2D mesh is already computed in ngMesh
1632 //================================================================================
1634 void NETGENPlugin_Mesher::AddIntVerticesInSolids(const netgen::OCCGeometry& occgeom,
1635 netgen::Mesh& ngMesh,
1636 vector<const SMDS_MeshNode*>& nodeVec,
1637 NETGENPlugin_Internals& internalShapes)
1639 #ifdef DUMP_TRIANGLES_SCRIPT
1640 // create a python script making a mesh containing triangles added for internal vertices
1641 ofstream py(DUMP_TRIANGLES_SCRIPT);
1642 py << "import SMESH"<< endl
1643 << "from salome.smesh import smeshBuilder"<<endl
1644 << "smesh = smeshBuilder.New(salome.myStudy)"<<endl
1645 << "m = smesh.Mesh(name='triangles')" << endl;
1647 if ((int) nodeVec.size() < ngMesh.GetNP() )
1648 nodeVec.resize( ngMesh.GetNP(), 0 );
1650 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1651 SMESH_MesherHelper helper( internalShapes.getMesh() );
1653 const map<int,list<int> >& so2Vert = internalShapes.getSolidsWithVertices();
1654 map<int,list<int> >::const_iterator s2v = so2Vert.begin();
1655 for ( ; s2v != so2Vert.end(); ++s2v )
1657 const TopoDS_Shape& solid = meshDS->IndexToShape( s2v->first );
1658 if ( solid.IsNull() ) continue;
1659 int solidNgID = occgeom.somap.FindIndex (solid);
1660 if ( solidNgID < 0 && !occgeom.somap.IsEmpty() ) continue;
1662 helper.SetSubShape( solid );
1663 helper.SetElementsOnShape( true );
1665 // find ng indices of faces within the solid
1667 for (TopExp_Explorer fExp(solid, TopAbs_FACE); fExp.More(); fExp.Next() )
1668 ngFaceIds.insert( occgeom.fmap.FindIndex( fExp.Current() ));
1669 if ( ngFaceIds.size() == 1 && *ngFaceIds.begin() == 0 )
1670 ngFaceIds.insert( 1 );
1672 // Get data of internal vertices and add them to ngMesh
1674 multimap< double, TIntVSoData > dist2VData; // sort vertices by distance from ng faces
1676 int i, nbFaceInit = ngMesh.GetNSE();
1678 // boundary characteristics
1679 double totSegLen = 0;
1682 const list<int>& iVertices = s2v->second;
1683 list<int>::const_iterator iv = iVertices.begin();
1684 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1687 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1689 // get node on vertex
1690 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1693 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1694 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1695 nV = SMESH_Algo::VertexNode( V, meshDS );
1696 if ( !nV ) continue;
1699 netgen::MeshPoint mpV( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1700 ngMesh.AddPoint ( mpV, 1, netgen::FIXEDPOINT );
1701 vData.ngId = ngMesh.GetNP();
1702 nodeVec.push_back( nV );
1704 // loop on all 2d elements to find the one closest to vertex and to count
1705 // average segment length
1706 double closeDist2 = numeric_limits<double>::max(), avgDist2;
1707 for (i = 1; i <= ngMesh.GetNSE(); ++i)
1709 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1710 if ( !ngFaceIds.count( elem.GetIndex() )) continue;
1712 multimap< double, int> dist2nID; // sort nodes of element by distance from V
1713 for ( int j = 0; j < elem.GetNP(); ++j)
1715 netgen::MeshPoint mp = ngMesh.Point( elem[j] );
1716 double d2 = dist2( mpV, mp );
1717 dist2nID.insert( make_pair( d2, elem[j] ));
1718 avgDist2 += d2 / elem.GetNP();
1720 totNbSeg++, totSegLen+= sqrt( dist2( mp, ngMesh.Point( elem[(j+1)%elem.GetNP()])));
1722 double dist = dist2nID.begin()->first; //avgDist2;
1723 if ( dist < closeDist2 )
1724 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= dist;
1726 dist2VData.insert( make_pair( closeDist2, vData ));
1729 if ( totNbSeg == 0 ) break;
1730 double avgSegLen = totSegLen / totNbSeg;
1732 // Loop on vertices to add triangles
1734 multimap< double, TIntVSoData >::iterator dist_vData = dist2VData.begin();
1735 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1737 double closeDist2 = dist_vData->first;
1738 TIntVSoData & vData = dist_vData->second;
1740 const netgen::MeshPoint& mpV = ngMesh.Point( vData.ngId );
1742 // try to find more close face among ones added for internal vertices
1743 for (i = nbFaceInit+1; i <= ngMesh.GetNSE(); ++i)
1745 double avgDist2 = 0;
1746 multimap< double, int> dist2nID;
1747 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1748 for ( int j = 0; j < elem.GetNP(); ++j)
1750 double d = dist2( mpV, ngMesh.Point( elem[j] ));
1751 dist2nID.insert( make_pair( d, elem[j] ));
1752 avgDist2 += d / elem.GetNP();
1753 if ( avgDist2 < closeDist2 )
1754 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= avgDist2;
1757 // sort nodes of the closest face by angle with vector from V to the closest node
1758 const double tol = numeric_limits<double>::min();
1759 map< double, int > angle2ID;
1760 const netgen::Element2d& closeFace = ngMesh.SurfaceElement( vData.ngIdClose );
1761 netgen::MeshPoint mp[2];
1762 mp[0] = ngMesh.Point( vData.ngIdCloseN );
1763 gp_XYZ p1( NGPOINT_COORDS( mp[0] ));
1764 gp_XYZ pV( NGPOINT_COORDS( mpV ));
1765 gp_Vec v2p1( pV, p1 );
1766 double distN1 = v2p1.Magnitude();
1767 if ( distN1 <= tol ) continue;
1769 for ( int j = 0; j < closeFace.GetNP(); ++j)
1771 mp[1] = ngMesh.Point( closeFace[j] );
1772 gp_Vec v2p( pV, gp_Pnt( NGPOINT_COORDS( mp[1] )) );
1773 angle2ID.insert( make_pair( v2p1.Angle( v2p ), closeFace[j]));
1775 // get node with angle of 60 degrees or greater
1776 map< double, int >::iterator angle_id = angle2ID.lower_bound( 60. * M_PI / 180. );
1777 if ( angle_id == angle2ID.end() ) angle_id = --angle2ID.end();
1778 const double minAngle = 30. * M_PI / 180.;
1779 const double angle = angle_id->first;
1780 bool angleOK = ( angle > minAngle );
1782 // find points to create a triangle
1783 netgen::Element2d tri(3);
1785 tri[0] = vData.ngId;
1786 tri[1] = vData.ngIdCloseN; // to use the closest nodes
1787 tri[2] = angle_id->second; // to use the node with best angle
1789 // decide whether to use the closest node and the node with best angle or to create new ones
1790 for ( int isBestAngleN = 0; isBestAngleN < 2; ++isBestAngleN )
1792 bool createNew = !angleOK; //, distOK = true;
1794 int triInd = isBestAngleN ? 2 : 1;
1795 mp[isBestAngleN] = ngMesh.Point( tri[triInd] );
1800 double distN2 = sqrt( dist2( mpV, mp[isBestAngleN]));
1801 createNew = ( fabs( distN2 - distN1 ) > 0.25 * distN1 );
1803 else if ( angle < tol )
1805 v2p1.SetX( v2p1.X() + 1e-3 );
1811 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1812 bool avgLenOK = ( avgSegLen < 0.75 * distN1 );
1813 bool hintLenOK = ( segLenHint < 0.75 * distN1 );
1814 createNew = (createNew || avgLenOK || hintLenOK );
1815 // we create a new node not closer than 0.5 to the closest face
1816 // in order not to clash with other close face
1817 double r = min( 0.5, ( hintLenOK ? segLenHint : avgSegLen ) / distN1 );
1818 distFromV = r * distN1;
1822 // create a new point, between the node and the vertex if angleOK
1823 gp_XYZ p( NGPOINT_COORDS( mp[isBestAngleN] ));
1824 gp_Vec v2p( pV, p ); v2p.Normalize();
1825 if ( isBestAngleN && !angleOK )
1826 p = p1 + gp_Dir( v2p.XYZ() - v2p1.XYZ()).XYZ() * distN1 * 0.95;
1828 p = pV + v2p.XYZ() * distFromV;
1830 if ( !isBestAngleN ) p1 = p, distN1 = distFromV;
1832 mp[isBestAngleN].SetPoint( netgen::Point<3> (p.X(), p.Y(), p.Z()));
1833 ngMesh.AddPoint ( mp[isBestAngleN], 1, netgen::SURFACEPOINT );
1834 tri[triInd] = ngMesh.GetNP();
1835 nodeVec.push_back( helper.AddNode( p.X(), p.Y(), p.Z()) );
1838 ngMesh.AddSurfaceElement (tri);
1839 swap( tri[1], tri[2] );
1840 ngMesh.AddSurfaceElement (tri);
1842 #ifdef DUMP_TRIANGLES_SCRIPT
1843 py << "n1 = m.AddNode( "<< mpV(0)<<", "<< mpV(1)<<", "<< mpV(2)<<") "<< endl
1844 << "n2 = m.AddNode( "<< mp[0](0)<<", "<< mp[0](1)<<", "<< mp[0](2)<<") "<< endl
1845 << "n3 = m.AddNode( "<< mp[1](0)<<", "<< mp[1](1)<<", "<< mp[1](2)<<" )" << endl
1846 << "m.AddFace([n1,n2,n3])" << endl;
1848 } // loop on internal vertices of a solid
1850 } // loop on solids with internal vertices
1853 //================================================================================
1855 * \brief Fill netgen mesh with segments of a FACE
1856 * \param ngMesh - netgen mesh
1857 * \param geom - container of OCCT geometry to mesh
1858 * \param wires - data of nodes on FACE boundary
1859 * \param helper - mesher helper holding the FACE
1860 * \param nodeVec - vector of nodes in which node index == netgen ID
1861 * \retval SMESH_ComputeErrorPtr - error description
1863 //================================================================================
1865 SMESH_ComputeErrorPtr
1866 NETGENPlugin_Mesher::AddSegmentsToMesh(netgen::Mesh& ngMesh,
1867 netgen::OCCGeometry& geom,
1868 const TSideVector& wires,
1869 SMESH_MesherHelper& helper,
1870 vector< const SMDS_MeshNode* > & nodeVec,
1871 const bool overrideMinH)
1873 // ----------------------------
1874 // Check wires and count nodes
1875 // ----------------------------
1877 for ( size_t iW = 0; iW < wires.size(); ++iW )
1879 StdMeshers_FaceSidePtr wire = wires[ iW ];
1880 if ( wire->MissVertexNode() )
1882 // Commented for issue 0020960. It worked for the case, let's wait for case where it doesn't.
1883 // It seems that there is no reason for this limitation
1885 // (new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH, "Missing nodes on vertices"));
1887 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1888 if ((int) uvPtVec.size() != wire->NbPoints() )
1889 return SMESH_ComputeError::New(COMPERR_BAD_INPUT_MESH,
1890 SMESH_Comment("Unexpected nb of points on wire ") << iW
1891 << ": " << uvPtVec.size()<<" != "<<wire->NbPoints());
1892 nbNodes += wire->NbPoints();
1894 nodeVec.reserve( nodeVec.size() + nbNodes + 1 );
1895 if ( nodeVec.empty() )
1896 nodeVec.push_back( 0 );
1898 // -----------------
1900 // -----------------
1902 const bool wasNgMeshEmpty = ( ngMesh.GetNP() < 1 ); /* true => this method is called by
1903 NETGENPlugin_NETGEN_2D_ONLY */
1905 // map for nodes on vertices since they can be shared between wires
1906 // ( issue 0020676, face_int_box.brep) and nodes built by NETGEN
1907 map<const SMDS_MeshNode*, int > node2ngID;
1908 if ( !wasNgMeshEmpty ) // fill node2ngID with nodes built by NETGEN
1910 set< int > subIDs; // ids of sub-shapes of the FACE
1911 for ( size_t iW = 0; iW < wires.size(); ++iW )
1913 StdMeshers_FaceSidePtr wire = wires[ iW ];
1914 for ( int iE = 0, nbE = wire->NbEdges(); iE < nbE; ++iE )
1916 subIDs.insert( wire->EdgeID( iE ));
1917 subIDs.insert( helper.GetMeshDS()->ShapeToIndex( wire->FirstVertex( iE )));
1920 for ( size_t ngID = 1; ngID < nodeVec.size(); ++ngID )
1921 if ( subIDs.count( nodeVec[ngID]->getshapeId() ))
1922 node2ngID.insert( make_pair( nodeVec[ngID], ngID ));
1925 const int solidID = 0, faceID = geom.fmap.FindIndex( helper.GetSubShape() );
1926 if ( ngMesh.GetNFD() < 1 )
1927 ngMesh.AddFaceDescriptor( netgen::FaceDescriptor( faceID, solidID, solidID, 0 ));
1929 for ( size_t iW = 0; iW < wires.size(); ++iW )
1931 StdMeshers_FaceSidePtr wire = wires[ iW ];
1932 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1933 const int nbSegments = wire->NbPoints() - 1;
1935 // assure the 1st node to be in node2ngID, which is needed to correctly
1936 // "close chain of segments" (see below) in case if the 1st node is not
1937 // onVertex because it is on a Viscous layer
1938 node2ngID.insert( make_pair( uvPtVec[ 0 ].node, ngMesh.GetNP() + 1 ));
1940 // compute length of every segment
1941 vector<double> segLen( nbSegments );
1942 for ( int i = 0; i < nbSegments; ++i )
1943 segLen[i] = SMESH_TNodeXYZ( uvPtVec[ i ].node ).Distance( uvPtVec[ i+1 ].node );
1945 int edgeID = 1, posID = -2;
1946 bool isInternalWire = false;
1947 double vertexNormPar = 0;
1948 const int prevNbNGSeg = ngMesh.GetNSeg();
1949 for ( int i = 0; i < nbSegments; ++i ) // loop on segments
1951 // Add the first point of a segment
1953 const SMDS_MeshNode * n = uvPtVec[ i ].node;
1954 const int posShapeID = n->getshapeId();
1955 bool onVertex = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX );
1956 bool onEdge = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_EDGE );
1958 // skip nodes on degenerated edges
1959 if ( helper.IsDegenShape( posShapeID ) &&
1960 helper.IsDegenShape( uvPtVec[ i+1 ].node->getshapeId() ))
1963 int ngID1 = ngMesh.GetNP() + 1, ngID2 = ngID1+1;
1964 if ( onVertex || ( !wasNgMeshEmpty && onEdge ) || helper.IsRealSeam( posShapeID ))
1965 ngID1 = node2ngID.insert( make_pair( n, ngID1 )).first->second;
1966 if ( ngID1 > ngMesh.GetNP() )
1968 netgen::MeshPoint mp( netgen::Point<3> (n->X(), n->Y(), n->Z()) );
1969 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1970 nodeVec.push_back( n );
1972 else // n is in ngMesh already, and ngID2 in prev segment is wrong
1974 ngID2 = ngMesh.GetNP() + 1;
1975 if ( i > 0 ) // prev segment belongs to same wire
1977 netgen::Segment& prevSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
1984 netgen::Segment seg;
1986 seg[0] = ngID1; // ng node id
1987 seg[1] = ngID2; // ng node id
1988 seg.edgenr = ngMesh.GetNSeg() + 1; // ng segment id
1989 seg.si = faceID; // = geom.fmap.FindIndex (face);
1991 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1993 const UVPtStruct& pnt = uvPtVec[ i + iEnd ];
1995 seg.epgeominfo[ iEnd ].dist = pnt.param; // param on curve
1996 seg.epgeominfo[ iEnd ].u = pnt.u;
1997 seg.epgeominfo[ iEnd ].v = pnt.v;
1999 // find out edge id and node parameter on edge
2000 onVertex = ( pnt.normParam + 1e-10 > vertexNormPar );
2001 if ( onVertex || posShapeID != posID )
2004 double normParam = pnt.normParam;
2006 normParam = 0.5 * ( uvPtVec[ i ].normParam + uvPtVec[ i+1 ].normParam );
2007 int edgeIndexInWire = wire->EdgeIndex( normParam );
2008 vertexNormPar = wire->LastParameter( edgeIndexInWire );
2009 const TopoDS_Edge& edge = wire->Edge( edgeIndexInWire );
2010 edgeID = geom.emap.FindIndex( edge );
2012 isInternalWire = ( edge.Orientation() == TopAbs_INTERNAL );
2013 // if ( onVertex ) // param on curve is different on each of two edges
2014 // seg.epgeominfo[ iEnd ].dist = helper.GetNodeU( edge, pnt.node );
2016 seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
2019 ngMesh.AddSegment (seg);
2021 // restrict size of elements near the segment
2022 SMESH_TNodeXYZ np1( n ), np2( uvPtVec[ i+1 ].node );
2023 // get an average size of adjacent segments to avoid sharp change of
2024 // element size (regression on issue 0020452, note 0010898)
2025 int iPrev = SMESH_MesherHelper::WrapIndex( i-1, nbSegments );
2026 int iNext = SMESH_MesherHelper::WrapIndex( i+1, nbSegments );
2027 double sumH = segLen[ iPrev ] + segLen[ i ] + segLen[ iNext ];
2028 int nbSeg = ( int( segLen[ iPrev ] > sumH / 100.) +
2029 int( segLen[ i ] > sumH / 100.) +
2030 int( segLen[ iNext ] > sumH / 100.));
2032 RestrictLocalSize( ngMesh, 0.5*(np1+np2), sumH / nbSeg, overrideMinH );
2034 if ( isInternalWire )
2036 swap (seg[0], seg[1]);
2037 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
2038 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
2039 ngMesh.AddSegment (seg);
2041 } // loop on segments on a wire
2043 // close chain of segments
2044 if ( nbSegments > 0 )
2046 netgen::Segment& lastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() - int( isInternalWire ));
2047 const SMDS_MeshNode * lastNode = uvPtVec.back().node;
2048 lastSeg[1] = node2ngID.insert( make_pair( lastNode, lastSeg[1] )).first->second;
2049 if ( lastSeg[1] > ngMesh.GetNP() )
2051 netgen::MeshPoint mp( netgen::Point<3> (lastNode->X(), lastNode->Y(), lastNode->Z()) );
2052 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
2053 nodeVec.push_back( lastNode );
2055 if ( isInternalWire )
2057 netgen::Segment& realLastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
2058 realLastSeg[0] = lastSeg[1];
2062 #ifdef DUMP_SEGMENTS
2063 cout << "BEGIN WIRE " << iW << endl;
2064 for ( int i = prevNbNGSeg+1; i <= ngMesh.GetNSeg(); ++i )
2066 netgen::Segment& seg = ngMesh.LineSegment( i );
2068 netgen::Segment& prevSeg = ngMesh.LineSegment( i-1 );
2069 if ( seg[0] == prevSeg[1] && seg[1] == prevSeg[0] )
2071 cout << "Segment: " << seg.edgenr << endl << "\tis REVRESE of the previous one" << endl;
2075 cout << "Segment: " << seg.edgenr << endl
2076 << "\tp1: " << seg[0] << " n" << nodeVec[ seg[0]]->GetID() << endl
2077 << "\tp2: " << seg[1] << " n" << nodeVec[ seg[1]]->GetID() << endl
2078 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
2079 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
2080 << "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
2081 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
2082 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl
2083 << "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
2085 cout << "--END WIRE " << iW << endl;
2087 SMESH_Comment __not_unused_variable( prevNbNGSeg );
2090 } // loop on WIREs of a FACE
2092 // add a segment instead of an internal vertex
2093 if ( wasNgMeshEmpty )
2095 NETGENPlugin_Internals intShapes( *helper.GetMesh(), helper.GetSubShape(), /*is3D=*/false );
2096 AddIntVerticesInFaces( geom, ngMesh, nodeVec, intShapes );
2098 ngMesh.CalcSurfacesOfNode();
2103 //================================================================================
2105 * \brief Fill SMESH mesh according to contents of netgen mesh
2106 * \param occgeo - container of OCCT geometry to mesh
2107 * \param ngMesh - netgen mesh
2108 * \param initState - bn of entities in netgen mesh before computing
2109 * \param sMesh - SMESH mesh to fill in
2110 * \param nodeVec - vector of nodes in which node index == netgen ID
2111 * \param comment - returns problem description
2112 * \param quadHelper - holder of medium nodes of sub-meshes
2113 * \retval int - error
2115 //================================================================================
2117 int NETGENPlugin_Mesher::FillSMesh(const netgen::OCCGeometry& occgeo,
2118 netgen::Mesh& ngMesh,
2119 const NETGENPlugin_ngMeshInfo& initState,
2121 std::vector<const SMDS_MeshNode*>& nodeVec,
2122 SMESH_Comment& comment,
2123 SMESH_MesherHelper* quadHelper)
2125 int nbNod = ngMesh.GetNP();
2126 int nbSeg = ngMesh.GetNSeg();
2127 int nbFac = ngMesh.GetNSE();
2128 int nbVol = ngMesh.GetNE();
2130 SMESHDS_Mesh* meshDS = sMesh.GetMeshDS();
2132 // quadHelper is used for either
2133 // 1) making quadratic elements when a lower dimention mesh is loaded
2134 // to SMESH before convertion to quadratic by NETGEN
2135 // 2) sewing of quadratic elements with quadratic elements of sub-meshes
2136 if ( quadHelper && !quadHelper->GetIsQuadratic() && quadHelper->GetTLinkNodeMap().empty() )
2139 // -------------------------------------
2140 // Create and insert nodes into nodeVec
2141 // -------------------------------------
2143 nodeVec.resize( nbNod + 1 );
2144 int i, nbInitNod = initState._nbNodes;
2145 for (i = nbInitNod+1; i <= nbNod; ++i )
2147 const netgen::MeshPoint& ngPoint = ngMesh.Point(i);
2148 SMDS_MeshNode* node = NULL;
2149 TopoDS_Vertex aVert;
2150 // First, netgen creates nodes on vertices in occgeo.vmap,
2151 // so node index corresponds to vertex index
2152 // but (issue 0020776) netgen does not create nodes with equal coordinates
2153 if ( i-nbInitNod <= occgeo.vmap.Extent() )
2155 gp_Pnt p ( NGPOINT_COORDS(ngPoint) );
2156 for (int iV = i-nbInitNod; aVert.IsNull() && iV <= occgeo.vmap.Extent(); ++iV)
2158 aVert = TopoDS::Vertex( occgeo.vmap( iV ));
2159 gp_Pnt pV = BRep_Tool::Pnt( aVert );
2160 if ( p.SquareDistance( pV ) > 1e-20 )
2163 node = const_cast<SMDS_MeshNode*>( SMESH_Algo::VertexNode( aVert, meshDS ));
2166 if (!node) // node not found on vertex
2168 node = meshDS->AddNode( NGPOINT_COORDS( ngPoint ));
2169 if (!aVert.IsNull())
2170 meshDS->SetNodeOnVertex(node, aVert);
2175 // -------------------------------------------
2176 // Create mesh segments along geometric edges
2177 // -------------------------------------------
2179 int nbInitSeg = initState._nbSegments;
2180 for (i = nbInitSeg+1; i <= nbSeg; ++i )
2182 const netgen::Segment& seg = ngMesh.LineSegment(i);
2184 int pinds[3] = { seg.pnums[0], seg.pnums[1], seg.pnums[2] };
2187 for (int j=0; j < 3; ++j)
2189 int pind = pinds[j];
2190 if (pind <= 0 || !nodeVec_ACCESS(pind))
2198 int aGeomEdgeInd = seg.epgeominfo[j].edgenr;
2199 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
2200 aEdge = TopoDS::Edge(occgeo.emap(aGeomEdgeInd));
2202 param = seg.epgeominfo[j].dist;
2205 else // middle point
2207 param = param2 * 0.5;
2209 if (!aEdge.IsNull() && nodeVec_ACCESS(pind)->getshapeId() < 1)
2211 meshDS->SetNodeOnEdge(nodeVec_ACCESS(pind), aEdge, param);
2216 SMDS_MeshEdge* edge = 0;
2217 if (nbp == 2) // second order ?
2219 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1])))
2221 if ( quadHelper ) // final mesh must be quadratic
2222 edge = quadHelper->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
2224 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
2228 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
2229 nodeVec_ACCESS(pinds[2])))
2231 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
2232 nodeVec_ACCESS(pinds[2]));
2236 if ( comment.empty() ) comment << "Cannot create a mesh edge";
2237 MESSAGE("Cannot create a mesh edge");
2238 nbSeg = nbFac = nbVol = 0;
2241 if ( !aEdge.IsNull() && edge->getshapeId() < 1 )
2242 meshDS->SetMeshElementOnShape(edge, aEdge);
2244 else if ( comment.empty() )
2246 comment << "Invalid netgen segment #" << i;
2250 // ----------------------------------------
2251 // Create mesh faces along geometric faces
2252 // ----------------------------------------
2254 int nbInitFac = initState._nbFaces;
2255 int quadFaceID = ngMesh.GetNFD() + 1;
2256 if ( nbInitFac < nbFac )
2257 // add a faces descriptor to exclude qudrangle elements generated by NETGEN
2258 // from computation of 3D mesh
2259 ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(quadFaceID, /*solid1=*/0, /*solid2=*/0, 0));
2261 vector<const SMDS_MeshNode*> nodes;
2262 for (i = nbInitFac+1; i <= nbFac; ++i )
2264 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
2265 const int aGeomFaceInd = elem.GetIndex();
2267 if (aGeomFaceInd > 0 && aGeomFaceInd <= occgeo.fmap.Extent())
2268 aFace = TopoDS::Face(occgeo.fmap(aGeomFaceInd));
2270 for ( int j = 1; j <= elem.GetNP(); ++j )
2272 int pind = elem.PNum(j);
2273 if ( pind < 1 || pind >= (int) nodeVec.size() )
2275 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind))
2277 nodes.push_back( node );
2278 if (!aFace.IsNull() && node->getshapeId() < 1)
2280 const netgen::PointGeomInfo& pgi = elem.GeomInfoPi(j);
2281 meshDS->SetNodeOnFace(node, aFace, pgi.u, pgi.v);
2285 if ((int) nodes.size() != elem.GetNP() )
2287 if ( comment.empty() )
2288 comment << "Invalid netgen 2d element #" << i;
2289 continue; // bad node ids
2291 SMDS_MeshFace* face = NULL;
2292 switch (elem.GetType())
2295 if ( quadHelper ) // final mesh must be quadratic
2296 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2]);
2298 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2]);
2301 if ( quadHelper ) // final mesh must be quadratic
2302 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2304 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2305 // exclude qudrangle elements from computation of 3D mesh
2306 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2309 nodes[5] = mediumNode( nodes[0],nodes[1],nodes[5], quadHelper );
2310 nodes[3] = mediumNode( nodes[1],nodes[2],nodes[3], quadHelper );
2311 nodes[4] = mediumNode( nodes[2],nodes[0],nodes[4], quadHelper );
2312 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[5],nodes[3],nodes[4]);
2315 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2316 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2317 nodes[5] = mediumNode( nodes[2],nodes[3],nodes[5], quadHelper );
2318 nodes[6] = mediumNode( nodes[3],nodes[0],nodes[6], quadHelper );
2319 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3],
2320 nodes[4],nodes[7],nodes[5],nodes[6]);
2321 // exclude qudrangle elements from computation of 3D mesh
2322 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2325 MESSAGE("NETGEN created a face of unexpected type, ignoring");
2330 if ( comment.empty() ) comment << "Cannot create a mesh face";
2331 MESSAGE("Cannot create a mesh face");
2332 nbSeg = nbFac = nbVol = 0;
2335 if ( !aFace.IsNull() )
2336 meshDS->SetMeshElementOnShape( face, aFace );
2339 // ------------------
2340 // Create tetrahedra
2341 // ------------------
2343 for ( i = 1; i <= nbVol; ++i )
2345 const netgen::Element& elem = ngMesh.VolumeElement(i);
2346 int aSolidInd = elem.GetIndex();
2347 TopoDS_Solid aSolid;
2348 if ( aSolidInd > 0 && aSolidInd <= occgeo.somap.Extent() )
2349 aSolid = TopoDS::Solid(occgeo.somap(aSolidInd));
2351 for ( int j = 1; j <= elem.GetNP(); ++j )
2353 int pind = elem.PNum(j);
2354 if ( pind < 1 || pind >= (int)nodeVec.size() )
2356 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind) )
2358 nodes.push_back(node);
2359 if ( !aSolid.IsNull() && node->getshapeId() < 1 )
2360 meshDS->SetNodeInVolume(node, aSolid);
2363 if ((int) nodes.size() != elem.GetNP() )
2365 if ( comment.empty() )
2366 comment << "Invalid netgen 3d element #" << i;
2369 SMDS_MeshVolume* vol = NULL;
2370 switch ( elem.GetType() )
2373 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3]);
2376 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2377 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2378 nodes[5] = mediumNode( nodes[2],nodes[0],nodes[5], quadHelper );
2379 nodes[6] = mediumNode( nodes[0],nodes[3],nodes[6], quadHelper );
2380 nodes[8] = mediumNode( nodes[1],nodes[3],nodes[8], quadHelper );
2381 nodes[9] = mediumNode( nodes[2],nodes[3],nodes[9], quadHelper );
2382 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3],
2383 nodes[4],nodes[7],nodes[5],nodes[6],nodes[8],nodes[9]);
2386 MESSAGE("NETGEN created a volume of unexpected type, ignoring");
2391 if ( comment.empty() ) comment << "Cannot create a mesh volume";
2392 MESSAGE("Cannot create a mesh volume");
2393 nbSeg = nbFac = nbVol = 0;
2396 if (!aSolid.IsNull())
2397 meshDS->SetMeshElementOnShape(vol, aSolid);
2399 return comment.empty() ? 0 : 1;
2404 //================================================================================
2406 * \brief Convert error into text
2408 //================================================================================
2410 std::string text(int err)
2415 SMESH_Comment("Error in netgen::OCCGenerateMesh() at ") << netgen::multithread.task;
2418 //================================================================================
2420 * \brief Convert exception into text
2422 //================================================================================
2424 std::string text(Standard_Failure& ex)
2426 SMESH_Comment str("Exception in netgen::OCCGenerateMesh()");
2427 str << " at " << netgen::multithread.task
2428 << ": " << ex.DynamicType()->Name();
2429 if ( ex.GetMessageString() && strlen( ex.GetMessageString() ))
2430 str << ": " << ex.GetMessageString();
2433 //================================================================================
2435 * \brief Convert exception into text
2437 //================================================================================
2439 std::string text(netgen::NgException& ex)
2441 SMESH_Comment str("NgException");
2442 if ( strlen( netgen::multithread.task ) > 0 )
2443 str << " at " << netgen::multithread.task;
2444 str << ": " << ex.What();
2448 //================================================================================
2450 * \brief Looks for triangles lying on a SOLID
2452 //================================================================================
2454 bool hasBadElemOnSolid( const list<const SMDS_MeshElement*>& elems,
2455 SMESH_subMesh* solidSM )
2457 TopTools_IndexedMapOfShape solidSubs;
2458 TopExp::MapShapes( solidSM->GetSubShape(), solidSubs );
2459 SMESHDS_Mesh* mesh = solidSM->GetFather()->GetMeshDS();
2461 list<const SMDS_MeshElement*>::const_iterator e = elems.begin();
2462 for ( ; e != elems.end(); ++e )
2464 const SMDS_MeshElement* elem = *e;
2465 // if ( elem->GetType() != SMDSAbs_Face ) -- 23047
2467 int nbNodesOnSolid = 0, nbNodes = elem->NbNodes();
2468 SMDS_NodeIteratorPtr nIt = elem->nodeIterator();
2469 while ( nIt->more() )
2471 const SMDS_MeshNode* n = nIt->next();
2472 const TopoDS_Shape& s = mesh->IndexToShape( n->getshapeId() );
2473 nbNodesOnSolid += ( !s.IsNull() && solidSubs.Contains( s ));
2474 if ( nbNodesOnSolid > 2 ||
2475 nbNodesOnSolid == nbNodes)
2482 const double edgeMeshingTime = 0.001;
2483 const double faceMeshingTime = 0.019;
2484 const double edgeFaceMeshingTime = edgeMeshingTime + faceMeshingTime;
2485 const double faceOptimizTime = 0.06;
2486 const double voluMeshingTime = 0.15;
2487 const double volOptimizeTime = 0.77;
2490 //=============================================================================
2492 * Here we are going to use the NETGEN mesher
2494 //=============================================================================
2496 bool NETGENPlugin_Mesher::Compute()
2498 NETGENPlugin_NetgenLibWrapper ngLib;
2500 netgen::MeshingParameters& mparams = netgen::mparam;
2502 SMESH_ComputeErrorPtr error = SMESH_ComputeError::New();
2503 SMESH_MesherHelper quadHelper( *_mesh );
2504 quadHelper.SetIsQuadratic( mparams.secondorder );
2506 // -------------------------
2507 // Prepare OCC geometry
2508 // -------------------------
2510 netgen::OCCGeometry occgeo;
2511 list< SMESH_subMesh* > meshedSM[3]; // for 0-2 dimensions
2512 NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
2513 PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
2516 _totalTime = edgeFaceMeshingTime;
2518 _totalTime += faceOptimizTime;
2520 _totalTime += voluMeshingTime + ( _optimize ? volOptimizeTime : 0 );
2521 double doneTime = 0;
2524 _curShapeIndex = -1;
2526 // -------------------------
2527 // Generate the mesh
2528 // -------------------------
2531 NETGENPlugin_ngMeshInfo initState; // it remembers size of ng mesh equal to size of Smesh
2533 SMESH_Comment comment;
2536 // vector of nodes in which node index == netgen ID
2537 vector< const SMDS_MeshNode* > nodeVec;
2545 // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
2546 mparams.uselocalh = false;
2547 mparams.grading = 0.8; // not limitited size growth
2549 if ( _simpleHyp->GetNumberOfSegments() )
2551 mparams.maxh = occgeo.boundingbox.Diam();
2554 mparams.maxh = _simpleHyp->GetLocalLength();
2557 if ( mparams.maxh == 0.0 )
2558 mparams.maxh = occgeo.boundingbox.Diam();
2559 if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
2560 mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
2562 // Local size on faces
2563 occgeo.face_maxh = mparams.maxh;
2565 // Let netgen create _ngMesh and calculate element size on not meshed shapes
2569 int startWith = netgen::MESHCONST_ANALYSE;
2570 int endWith = netgen::MESHCONST_ANALYSE;
2575 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2577 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2579 if(netgen::multithread.terminate)
2582 comment << text(err);
2584 catch (Standard_Failure& ex)
2586 comment << text(ex);
2588 catch (netgen::NgException & ex)
2590 comment << text(ex);
2591 if ( mparams.meshsizefilename )
2592 throw SMESH_ComputeError(COMPERR_BAD_PARMETERS, comment );
2594 err = 0; //- MESHCONST_ANALYSE isn't so important step
2597 ngLib.setMesh(( Ng_Mesh*) _ngMesh );
2599 _ngMesh->ClearFaceDescriptors(); // we make descriptors our-self
2601 if ( !mparams.uselocalh ) // mparams.grading is not taken into account yet
2602 _ngMesh->LocalHFunction().SetGrading( mparams.grading );
2606 // Pass 1D simple parameters to NETGEN
2607 // --------------------------------
2608 int nbSeg = _simpleHyp->GetNumberOfSegments();
2609 double segSize = _simpleHyp->GetLocalLength();
2610 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
2612 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
2614 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
2615 setLocalSize( e, segSize, *_ngMesh );
2618 else // if ( ! _simpleHyp )
2620 // Local size on shapes
2621 SetLocalSize( occgeo, *_ngMesh );
2624 // Precompute internal edges (issue 0020676) in order to
2625 // add mesh on them correctly (twice) to netgen mesh
2626 if ( !err && internals.hasInternalEdges() )
2628 // load internal shapes into OCCGeometry
2629 netgen::OCCGeometry intOccgeo;
2630 internals.getInternalEdges( intOccgeo.fmap, intOccgeo.emap, intOccgeo.vmap, meshedSM );
2631 intOccgeo.boundingbox = occgeo.boundingbox;
2632 intOccgeo.shape = occgeo.shape;
2633 intOccgeo.face_maxh.SetSize(intOccgeo.fmap.Extent());
2634 intOccgeo.face_maxh = netgen::mparam.maxh;
2635 netgen::Mesh *tmpNgMesh = NULL;
2639 // compute local H on internal shapes in the main mesh
2640 //OCCSetLocalMeshSize(intOccgeo, *_ngMesh); it deletes _ngMesh->localH
2642 // let netgen create a temporary mesh
2644 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2646 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2648 if(netgen::multithread.terminate)
2651 // copy LocalH from the main to temporary mesh
2652 initState.transferLocalH( _ngMesh, tmpNgMesh );
2654 // compute mesh on internal edges
2655 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2657 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2659 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2661 comment << text(err);
2663 catch (Standard_Failure& ex)
2665 comment << text(ex);
2668 initState.restoreLocalH( tmpNgMesh );
2670 // fill SMESH by netgen mesh
2671 vector< const SMDS_MeshNode* > tmpNodeVec;
2672 FillSMesh( intOccgeo, *tmpNgMesh, initState, *_mesh, tmpNodeVec, comment );
2673 err = ( err || !comment.empty() );
2675 nglib::Ng_DeleteMesh((nglib::Ng_Mesh*)tmpNgMesh);
2678 // Fill _ngMesh with nodes and segments of computed submeshes
2681 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_0D ]) &&
2682 FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_1D ], &quadHelper));
2684 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2689 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2694 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2696 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2698 if(netgen::multithread.terminate)
2701 comment << text(err);
2703 catch (Standard_Failure& ex)
2705 comment << text(ex);
2710 _ticTime = ( doneTime += edgeMeshingTime ) / _totalTime / _progressTic;
2712 mparams.uselocalh = true; // restore as it is used at surface optimization
2714 // ---------------------
2715 // compute surface mesh
2716 // ---------------------
2719 // Pass 2D simple parameters to NETGEN
2721 if ( double area = _simpleHyp->GetMaxElementArea() ) {
2723 mparams.maxh = sqrt(2. * area/sqrt(3.0));
2724 mparams.grading = 0.4; // moderate size growth
2727 // length from edges
2728 if ( _ngMesh->GetNSeg() ) {
2729 double edgeLength = 0;
2730 TopTools_MapOfShape visitedEdges;
2731 for ( TopExp_Explorer exp( _shape, TopAbs_EDGE ); exp.More(); exp.Next() )
2732 if( visitedEdges.Add(exp.Current()) )
2733 edgeLength += SMESH_Algo::EdgeLength( TopoDS::Edge( exp.Current() ));
2734 // we have to multiply length by 2 since for each TopoDS_Edge there
2735 // are double set of NETGEN edges, in other words, we have to
2736 // divide _ngMesh->GetNSeg() by 2.
2737 mparams.maxh = 2*edgeLength / _ngMesh->GetNSeg();
2740 mparams.maxh = 1000;
2742 mparams.grading = 0.2; // slow size growth
2744 mparams.quad = _simpleHyp->GetAllowQuadrangles();
2745 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2746 _ngMesh->SetGlobalH (mparams.maxh);
2747 netgen::Box<3> bb = occgeo.GetBoundingBox();
2748 bb.Increase (bb.Diam()/20);
2749 _ngMesh->SetLocalH (bb.PMin(), bb.PMax(), mparams.grading);
2752 // Care of vertices internal in faces (issue 0020676)
2753 if ( internals.hasInternalVertexInFace() )
2755 // store computed segments in SMESH in order not to create SMESH
2756 // edges for ng segments added by AddIntVerticesInFaces()
2757 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2758 // add segments to faces with internal vertices
2759 AddIntVerticesInFaces( occgeo, *_ngMesh, nodeVec, internals );
2760 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2763 // Build viscous layers
2764 if ( _isViscousLayers2D ||
2765 StdMeshers_ViscousLayers2D::HasProxyMesh( TopoDS::Face( occgeo.fmap(1) ), *_mesh ))
2767 if ( !internals.hasInternalVertexInFace() ) {
2768 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2769 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2771 SMESH_ProxyMesh::Ptr viscousMesh;
2772 SMESH_MesherHelper helper( *_mesh );
2773 for ( int faceID = 1; faceID <= occgeo.fmap.Extent(); ++faceID )
2775 const TopoDS_Face& F = TopoDS::Face( occgeo.fmap( faceID ));
2776 viscousMesh = StdMeshers_ViscousLayers2D::Compute( *_mesh, F );
2779 if ( viscousMesh->NbProxySubMeshes() == 0 )
2781 // exclude from computation ng segments built on EDGEs of F
2782 for (int i = 1; i <= _ngMesh->GetNSeg(); i++)
2784 netgen::Segment & seg = _ngMesh->LineSegment(i);
2785 if (seg.si == faceID)
2788 // add new segments to _ngMesh instead of excluded ones
2789 helper.SetSubShape( F );
2791 StdMeshers_FaceSide::GetFaceWires( F, *_mesh, /*skipMediumNodes=*/true,
2792 error, viscousMesh );
2793 error = AddSegmentsToMesh( *_ngMesh, occgeo, wires, helper, nodeVec );
2795 if ( !error ) error = SMESH_ComputeError::New();
2797 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2800 // Let netgen compute 2D mesh
2801 startWith = netgen::MESHCONST_MESHSURFACE;
2802 endWith = _optimize ? netgen::MESHCONST_OPTSURFACE : netgen::MESHCONST_MESHSURFACE;
2807 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2809 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2811 if(netgen::multithread.terminate)
2814 comment << text (err);
2816 catch (Standard_Failure& ex)
2818 comment << text(ex);
2819 //err = 1; -- try to make volumes anyway
2821 catch (netgen::NgException exc)
2823 comment << text(exc);
2824 //err = 1; -- try to make volumes anyway
2829 doneTime += faceMeshingTime + ( _optimize ? faceOptimizTime : 0 );
2830 _ticTime = doneTime / _totalTime / _progressTic;
2832 // ---------------------
2833 // generate volume mesh
2834 // ---------------------
2835 // Fill _ngMesh with nodes and faces of computed 2D submeshes
2836 if ( !err && _isVolume && ( !meshedSM[ MeshDim_2D ].empty() || mparams.quad ))
2838 // load SMESH with computed segments and faces
2839 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2841 // compute pyramids on quadrangles
2842 SMESH_ProxyMesh::Ptr proxyMesh;
2843 if ( _mesh->NbQuadrangles() > 0 )
2844 for ( int iS = 1; iS <= occgeo.somap.Extent(); ++iS )
2846 StdMeshers_QuadToTriaAdaptor* Adaptor = new StdMeshers_QuadToTriaAdaptor;
2847 proxyMesh.reset( Adaptor );
2849 int nbPyrams = _mesh->NbPyramids();
2850 Adaptor->Compute( *_mesh, occgeo.somap(iS) );
2851 if ( nbPyrams != _mesh->NbPyramids() )
2853 list< SMESH_subMesh* > quadFaceSM;
2854 for (TopExp_Explorer face(occgeo.somap(iS), TopAbs_FACE); face.More(); face.Next())
2855 if ( Adaptor->GetProxySubMesh( face.Current() ))
2857 quadFaceSM.push_back( _mesh->GetSubMesh( face.Current() ));
2858 meshedSM[ MeshDim_2D ].remove( quadFaceSM.back() );
2860 FillNgMesh(occgeo, *_ngMesh, nodeVec, quadFaceSM, &quadHelper, proxyMesh);
2863 // fill _ngMesh with faces of sub-meshes
2864 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_2D ], &quadHelper));
2865 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2866 //toPython( _ngMesh, "/tmp/ngPython.py");
2868 if (!err && _isVolume)
2870 // Pass 3D simple parameters to NETGEN
2871 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
2872 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
2874 if ( double vol = simple3d->GetMaxElementVolume() ) {
2876 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
2877 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2880 // length from faces
2881 mparams.maxh = _ngMesh->AverageH();
2883 _ngMesh->SetGlobalH (mparams.maxh);
2884 mparams.grading = 0.4;
2886 _ngMesh->CalcLocalH(mparams.grading);
2888 _ngMesh->CalcLocalH();
2891 // Care of vertices internal in solids and internal faces (issue 0020676)
2892 if ( internals.hasInternalVertexInSolid() || internals.hasInternalFaces() )
2894 // store computed faces in SMESH in order not to create SMESH
2895 // faces for ng faces added here
2896 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2897 // add ng faces to solids with internal vertices
2898 AddIntVerticesInSolids( occgeo, *_ngMesh, nodeVec, internals );
2899 // duplicate mesh faces on internal faces
2900 FixIntFaces( occgeo, *_ngMesh, internals );
2901 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2903 // Let netgen compute 3D mesh
2904 startWith = endWith = netgen::MESHCONST_MESHVOLUME;
2909 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2911 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2913 if(netgen::multithread.terminate)
2916 if ( comment.empty() ) // do not overwrite a previos error
2917 comment << text(err);
2919 catch (Standard_Failure& ex)
2921 if ( comment.empty() ) // do not overwrite a previos error
2922 comment << text(ex);
2925 catch (netgen::NgException exc)
2927 if ( comment.empty() ) // do not overwrite a previos error
2928 comment << text(exc);
2931 _ticTime = ( doneTime += voluMeshingTime ) / _totalTime / _progressTic;
2933 // Let netgen optimize 3D mesh
2934 if ( !err && _optimize )
2936 startWith = endWith = netgen::MESHCONST_OPTVOLUME;
2941 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2943 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2945 if(netgen::multithread.terminate)
2948 if ( comment.empty() ) // do not overwrite a previos error
2949 comment << text(err);
2951 catch (Standard_Failure& ex)
2953 if ( comment.empty() ) // do not overwrite a previos error
2954 comment << text(ex);
2956 catch (netgen::NgException exc)
2958 if ( comment.empty() ) // do not overwrite a previos error
2959 comment << text(exc);
2963 if (!err && mparams.secondorder > 0)
2968 if ( !meshedSM[ MeshDim_1D ].empty() )
2970 // remove segments not attached to geometry (IPAL0052479)
2971 for (int i = 1; i <= _ngMesh->GetNSeg(); ++i)
2973 const netgen::Segment & seg = _ngMesh->LineSegment (i);
2974 if ( seg.epgeominfo[ 0 ].edgenr == 0 )
2975 _ngMesh->DeleteSegment( i );
2977 _ngMesh->Compress();
2979 // convert to quadratic
2980 netgen::OCCRefinementSurfaces ref (occgeo);
2981 ref.MakeSecondOrder (*_ngMesh);
2983 // care of elements already loaded to SMESH
2984 // if ( initState._nbSegments > 0 )
2985 // makeQuadratic( occgeo.emap, _mesh );
2986 // if ( initState._nbFaces > 0 )
2987 // makeQuadratic( occgeo.fmap, _mesh );
2989 catch (Standard_Failure& ex)
2991 if ( comment.empty() ) // do not overwrite a previos error
2992 comment << "Exception in netgen at passing to 2nd order ";
2994 catch (netgen::NgException exc)
2996 if ( comment.empty() ) // do not overwrite a previos error
2997 comment << exc.What();
3002 _ticTime = 0.98 / _progressTic;
3004 //int nbNod = _ngMesh->GetNP();
3005 //int nbSeg = _ngMesh->GetNSeg();
3006 int nbFac = _ngMesh->GetNSE();
3007 int nbVol = _ngMesh->GetNE();
3008 bool isOK = ( !err && (_isVolume ? (nbVol > 0) : (nbFac > 0)) );
3010 // Feed back the SMESHDS with the generated Nodes and Elements
3011 if ( true /*isOK*/ ) // get whatever built
3013 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
3015 if ( quadHelper.GetIsQuadratic() ) // remove free nodes
3016 for ( size_t i = 0; i < nodeVec.size(); ++i )
3017 if ( nodeVec[i] && nodeVec[i]->NbInverseElements() == 0 )
3018 _mesh->GetMeshDS()->RemoveFreeNode( nodeVec[i], 0, /*fromGroups=*/false );
3020 SMESH_ComputeErrorPtr readErr = ReadErrors(nodeVec);
3021 if ( readErr && !readErr->myBadElements.empty() )
3024 if ( !comment.empty() && !readErr->myComment.empty() ) comment += "\n";
3025 comment += readErr->myComment;
3027 if ( error->IsOK() && ( !isOK || comment.size() > 0 ))
3028 error->myName = COMPERR_ALGO_FAILED;
3029 if ( !comment.empty() )
3030 error->myComment = comment;
3032 // SetIsAlwaysComputed( true ) to empty sub-meshes, which
3033 // appear if the geometry contains coincident sub-shape due
3034 // to bool merge_solids = 1; in netgen/libsrc/occ/occgenmesh.cpp
3035 const int nbMaps = 2;
3036 const TopTools_IndexedMapOfShape* geoMaps[nbMaps] =
3037 { & occgeo.vmap, & occgeo.emap/*, & occgeo.fmap*/ };
3038 for ( int iMap = 0; iMap < nbMaps; ++iMap )
3039 for (int i = 1; i <= geoMaps[iMap]->Extent(); i++)
3040 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( geoMaps[iMap]->FindKey(i)))
3041 if ( !sm->IsMeshComputed() )
3042 sm->SetIsAlwaysComputed( true );
3044 // set bad compute error to subshapes of all failed sub-shapes
3045 if ( !error->IsOK() )
3047 bool pb2D = false, pb3D = false;
3048 for (int i = 1; i <= occgeo.fmap.Extent(); i++) {
3049 int status = occgeo.facemeshstatus[i-1];
3050 if (status == netgen::FACE_MESHED_OK ) continue;
3051 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.fmap( i ))) {
3052 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
3053 if ( !smError || smError->IsOK() ) {
3054 if ( status == netgen::FACE_FAILED )
3055 smError.reset( new SMESH_ComputeError( *error ));
3057 smError.reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED, "Ignored" ));
3058 if ( SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
3059 smError->myName = COMPERR_WARNING;
3061 pb2D = pb2D || smError->IsKO();
3064 if ( !pb2D ) // all faces are OK
3065 for (int i = 1; i <= occgeo.somap.Extent(); i++)
3066 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.somap( i )))
3068 bool smComputed = nbVol && !sm->IsEmpty();
3069 if ( smComputed && internals.hasInternalVertexInSolid( sm->GetId() ))
3071 int nbIntV = internals.getSolidsWithVertices().find( sm->GetId() )->second.size();
3072 SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
3073 smComputed = ( smDS->NbElements() > 0 || smDS->NbNodes() > nbIntV );
3075 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
3076 if ( !smComputed && ( !smError || smError->IsOK() ))
3078 smError.reset( new SMESH_ComputeError( *error ));
3079 if ( nbVol && SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
3081 smError->myName = COMPERR_WARNING;
3083 else if ( !smError->myBadElements.empty() ) // bad surface mesh
3085 if ( !hasBadElemOnSolid( smError->myBadElements, sm ))
3089 pb3D = pb3D || ( smError && smError->IsKO() );
3091 if ( !pb2D && !pb3D )
3092 err = 0; // no fatal errors, only warnings
3095 ngLib._isComputeOk = !err;
3100 //=============================================================================
3104 //=============================================================================
3105 bool NETGENPlugin_Mesher::Evaluate(MapShapeNbElems& aResMap)
3107 netgen::MeshingParameters& mparams = netgen::mparam;
3110 // -------------------------
3111 // Prepare OCC geometry
3112 // -------------------------
3113 netgen::OCCGeometry occgeo;
3114 list< SMESH_subMesh* > meshedSM[4]; // for 0-3 dimensions
3115 NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
3116 PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
3118 bool tooManyElems = false;
3119 const int hugeNb = std::numeric_limits<int>::max() / 100;
3124 // pass 1D simple parameters to NETGEN
3127 // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
3128 mparams.uselocalh = false;
3129 mparams.grading = 0.8; // not limitited size growth
3131 if ( _simpleHyp->GetNumberOfSegments() )
3133 mparams.maxh = occgeo.boundingbox.Diam();
3136 mparams.maxh = _simpleHyp->GetLocalLength();
3139 if ( mparams.maxh == 0.0 )
3140 mparams.maxh = occgeo.boundingbox.Diam();
3141 if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
3142 mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
3144 // let netgen create _ngMesh and calculate element size on not meshed shapes
3145 NETGENPlugin_NetgenLibWrapper ngLib;
3146 netgen::Mesh *ngMesh = NULL;
3150 int startWith = netgen::MESHCONST_ANALYSE;
3151 int endWith = netgen::MESHCONST_MESHEDGES;
3153 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, mparams, startWith, endWith);
3155 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, startWith, endWith, optstr);
3158 if(netgen::multithread.terminate)
3161 ngLib.setMesh(( Ng_Mesh*) ngMesh );
3163 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( _shape ))
3164 sm->GetComputeError().reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED ));
3169 // Pass 1D simple parameters to NETGEN
3170 // --------------------------------
3171 int nbSeg = _simpleHyp->GetNumberOfSegments();
3172 double segSize = _simpleHyp->GetLocalLength();
3173 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
3175 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
3177 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
3178 setLocalSize( e, segSize, *ngMesh );
3181 else // if ( ! _simpleHyp )
3183 // Local size on shapes
3184 SetLocalSize( occgeo, *ngMesh );
3186 // calculate total nb of segments and length of edges
3187 double fullLen = 0.0;
3189 int entity = mparams.secondorder > 0 ? SMDSEntity_Quad_Edge : SMDSEntity_Edge;
3190 TopTools_DataMapOfShapeInteger Edge2NbSeg;
3191 for (TopExp_Explorer exp(_shape, TopAbs_EDGE); exp.More(); exp.Next())
3193 TopoDS_Edge E = TopoDS::Edge( exp.Current() );
3194 if( !Edge2NbSeg.Bind(E,0) )
3197 double aLen = SMESH_Algo::EdgeLength(E);
3200 vector<int>& aVec = aResMap[_mesh->GetSubMesh(E)];
3202 aVec.resize( SMDSEntity_Last, 0);
3204 fullNbSeg += aVec[ entity ];
3207 // store nb of segments computed by Netgen
3208 NCollection_Map<Link> linkMap;
3209 for (int i = 1; i <= ngMesh->GetNSeg(); ++i )
3211 const netgen::Segment& seg = ngMesh->LineSegment(i);
3212 Link link(seg[0], seg[1]);
3213 if ( !linkMap.Add( link )) continue;
3214 int aGeomEdgeInd = seg.epgeominfo[0].edgenr;
3215 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
3217 vector<int>& aVec = aResMap[_mesh->GetSubMesh(occgeo.emap(aGeomEdgeInd))];
3221 // store nb of nodes on edges computed by Netgen
3222 TopTools_DataMapIteratorOfDataMapOfShapeInteger Edge2NbSegIt(Edge2NbSeg);
3223 for (; Edge2NbSegIt.More(); Edge2NbSegIt.Next())
3225 vector<int>& aVec = aResMap[_mesh->GetSubMesh(Edge2NbSegIt.Key())];
3226 if ( aVec[ entity ] > 1 && aVec[ SMDSEntity_Node ] == 0 )
3227 aVec[SMDSEntity_Node] = mparams.secondorder > 0 ? 2*aVec[ entity ]-1 : aVec[ entity ]-1;
3229 fullNbSeg += aVec[ entity ];
3230 Edge2NbSeg( Edge2NbSegIt.Key() ) = aVec[ entity ];
3232 if ( fullNbSeg == 0 )
3239 if ( double area = _simpleHyp->GetMaxElementArea() ) {
3241 mparams.maxh = sqrt(2. * area/sqrt(3.0));
3242 mparams.grading = 0.4; // moderate size growth
3245 // length from edges
3246 mparams.maxh = fullLen/fullNbSeg;
3247 mparams.grading = 0.2; // slow size growth
3250 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3251 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3253 for (TopExp_Explorer exp(_shape, TopAbs_FACE); exp.More(); exp.Next())
3255 TopoDS_Face F = TopoDS::Face( exp.Current() );
3256 SMESH_subMesh *sm = _mesh->GetSubMesh(F);
3258 BRepGProp::SurfaceProperties(F,G);
3259 double anArea = G.Mass();
3260 tooManyElems = tooManyElems || ( anArea/hugeNb > mparams.maxh*mparams.maxh );
3262 if ( !tooManyElems )
3264 TopTools_MapOfShape egdes;
3265 for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next())
3266 if ( egdes.Add( exp1.Current() ))
3267 nb1d += Edge2NbSeg.Find(exp1.Current());
3269 int nbFaces = tooManyElems ? hugeNb : int( 4*anArea / (mparams.maxh*mparams.maxh*sqrt(3.)));
3270 int nbNodes = tooManyElems ? hugeNb : (( nbFaces*3 - (nb1d-1)*2 ) / 6 + 1 );
3272 vector<int> aVec(SMDSEntity_Last, 0);
3273 if( mparams.secondorder > 0 ) {
3274 int nb1d_in = (nbFaces*3 - nb1d) / 2;
3275 aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
3276 aVec[SMDSEntity_Quad_Triangle] = nbFaces;
3279 aVec[SMDSEntity_Node] = Max ( nbNodes, 0 );
3280 aVec[SMDSEntity_Triangle] = nbFaces;
3282 aResMap[sm].swap(aVec);
3289 // pass 3D simple parameters to NETGEN
3290 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
3291 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
3293 if ( double vol = simple3d->GetMaxElementVolume() ) {
3295 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
3296 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3299 // using previous length from faces
3301 mparams.grading = 0.4;
3302 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3305 BRepGProp::VolumeProperties(_shape,G);
3306 double aVolume = G.Mass();
3307 double tetrVol = 0.1179*mparams.maxh*mparams.maxh*mparams.maxh;
3308 tooManyElems = tooManyElems || ( aVolume/hugeNb > tetrVol );
3309 int nbVols = tooManyElems ? hugeNb : int(aVolume/tetrVol);
3310 int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
3311 vector<int> aVec(SMDSEntity_Last, 0 );
3312 if ( tooManyElems ) // avoid FPE
3314 aVec[SMDSEntity_Node] = hugeNb;
3315 aVec[ mparams.secondorder > 0 ? SMDSEntity_Quad_Tetra : SMDSEntity_Tetra] = hugeNb;
3319 if( mparams.secondorder > 0 ) {
3320 aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
3321 aVec[SMDSEntity_Quad_Tetra] = nbVols;
3324 aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
3325 aVec[SMDSEntity_Tetra] = nbVols;
3328 SMESH_subMesh *sm = _mesh->GetSubMesh(_shape);
3329 aResMap[sm].swap(aVec);
3335 double NETGENPlugin_Mesher::GetProgress(const SMESH_Algo* holder,
3336 const int * algoProgressTic,
3337 const double * algoProgress) const
3339 ((int&) _progressTic ) = *algoProgressTic + 1;
3341 if ( !_occgeom ) return 0;
3343 double progress = -1;
3346 if ( _ticTime < 0 && netgen::multithread.task[0] == 'O'/*Optimizing surface*/ )
3348 ((double&) _ticTime ) = edgeFaceMeshingTime / _totalTime / _progressTic;
3350 else if ( !_optimize /*&& _occgeom->fmap.Extent() > 1*/ )
3352 int doneShapeIndex = -1;
3353 while ( doneShapeIndex+1 < _occgeom->facemeshstatus.Size() &&
3354 _occgeom->facemeshstatus[ doneShapeIndex+1 ])
3356 if ( doneShapeIndex+1 != _curShapeIndex )
3358 ((int&) _curShapeIndex) = doneShapeIndex+1;
3359 double doneShapeRate = _curShapeIndex / double( _occgeom->fmap.Extent() );
3360 double doneTime = edgeMeshingTime + doneShapeRate * faceMeshingTime;
3361 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3362 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3363 // << " " << doneTime / _totalTime / _progressTic << endl;
3367 else if ( !_optimize && _occgeom->somap.Extent() > 1 )
3369 int curShapeIndex = _curShapeIndex;
3370 if ( _ngMesh->GetNE() > 0 )
3372 netgen::Element el = (*_ngMesh)[netgen::ElementIndex( _ngMesh->GetNE()-1 )];
3373 curShapeIndex = el.GetIndex();
3375 if ( curShapeIndex != _curShapeIndex )
3377 ((int&) _curShapeIndex) = curShapeIndex;
3378 double doneShapeRate = _curShapeIndex / double( _occgeom->somap.Extent() );
3379 double doneTime = edgeFaceMeshingTime + doneShapeRate * voluMeshingTime;
3380 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3381 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3382 // << " " << doneTime / _totalTime / _progressTic << endl;
3387 progress = Max( *algoProgressTic * _ticTime, *algoProgress );
3392 netgen::multithread.task[0] == 'D'/*elaunay meshing*/ &&
3393 progress > voluMeshingTime )
3395 progress = voluMeshingTime;
3396 ((double&) _ticTime) = voluMeshingTime / _totalTime / _progressTic;
3398 ((int&) *algoProgressTic )++;
3399 ((double&) *algoProgress) = progress;
3401 //cout << progress << " " << *algoProgressTic << " " << netgen::multithread.task << " "<< _ticTime << endl;
3403 return Min( progress, 0.99 );
3406 //================================================================================
3408 * \brief Read mesh entities preventing successful computation from "test.out" file
3410 //================================================================================
3412 SMESH_ComputeErrorPtr
3413 NETGENPlugin_Mesher::ReadErrors(const vector<const SMDS_MeshNode* >& nodeVec)
3415 SMESH_ComputeErrorPtr err = SMESH_ComputeError::New
3416 (COMPERR_BAD_INPUT_MESH, "Some edges multiple times in surface mesh");
3417 SMESH_File file("test.out");
3419 vector<int> three1(3), three2(3);
3420 const char* badEdgeStr = " multiple times in surface mesh";
3421 const int badEdgeStrLen = strlen( badEdgeStr );
3422 const int nbNodes = nodeVec.size();
3424 while( !file.eof() )
3426 if ( strncmp( file, "Edge ", 5 ) == 0 &&
3427 file.getInts( two ) &&
3428 strncmp( file, badEdgeStr, badEdgeStrLen ) == 0 &&
3429 two[0] < nbNodes && two[1] < nbNodes )
3431 err->myBadElements.push_back( new SMDS_LinearEdge( nodeVec[ two[0]], nodeVec[ two[1]] ));
3432 file += badEdgeStrLen;
3434 else if ( strncmp( file, "Intersecting: ", 14 ) == 0 )
3437 // openelement 18 with open element 126
3441 const char* pos = file;
3442 bool ok = ( strncmp( file, "openelement ", 12 ) == 0 );
3443 ok = ok && file.getInts( two );
3444 ok = ok && file.getInts( three1 );
3445 ok = ok && file.getInts( three2 );
3446 for ( int i = 0; ok && i < 3; ++i )
3447 ok = ( three1[i] < nbNodes && nodeVec[ three1[i]]);
3448 for ( int i = 0; ok && i < 3; ++i )
3449 ok = ( three2[i] < nbNodes && nodeVec[ three2[i]]);
3452 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three1[0]],
3453 nodeVec[ three1[1]],
3454 nodeVec[ three1[2]]));
3455 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three2[0]],
3456 nodeVec[ three2[1]],
3457 nodeVec[ three2[2]]));
3458 err->myComment = "Intersecting triangles";
3472 size_t nbBadElems = err->myBadElements.size();
3473 if ( nbBadElems ) nbBadElems++; // avoid warning: variable set but not used
3479 //================================================================================
3481 * \brief Write a python script creating an equivalent SALOME mesh.
3482 * This is useful to see what mesh is passed as input for the next step of mesh
3483 * generation (of mesh of higher dimension)
3485 //================================================================================
3487 void NETGENPlugin_Mesher::toPython( const netgen::Mesh* ngMesh )
3489 const char* pyFile = "/tmp/ngMesh.py";
3490 ofstream outfile( pyFile, ios::out );
3491 if ( !outfile ) return;
3493 outfile << "import salome, SMESH" << endl
3494 << "from salome.smesh import smeshBuilder" << endl
3495 << "smesh = smeshBuilder.New(salome.myStudy)" << endl
3496 << "mesh = smesh.Mesh()" << endl << endl;
3498 using namespace netgen;
3500 for (pi = PointIndex::BASE;
3501 pi < ngMesh->GetNP()+PointIndex::BASE; pi++)
3503 outfile << "mesh.AddNode( ";
3504 outfile << (*ngMesh)[pi](0) << ", ";
3505 outfile << (*ngMesh)[pi](1) << ", ";
3506 outfile << (*ngMesh)[pi](2) << ") ## "<< pi << endl;
3509 int nbDom = ngMesh->GetNDomains();
3510 for ( int i = 0; i < nbDom; ++i )
3511 outfile<< "grp" << i+1 << " = mesh.CreateEmptyGroup( SMESH.FACE, 'domain"<< i+1 << "')"<< endl;
3513 SurfaceElementIndex sei;
3514 for (sei = 0; sei < ngMesh->GetNSE(); sei++)
3516 outfile << "mesh.AddFace([ ";
3517 Element2d sel = (*ngMesh)[sei];
3518 for (int j = 0; j < sel.GetNP(); j++)
3519 outfile << sel[j] << ( j+1 < sel.GetNP() ? ", " : " ])");
3520 if ( sel.IsDeleted() ) outfile << " ## IsDeleted ";
3523 if ((*ngMesh)[sei].GetIndex())
3525 if ( int dom1 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainIn())
3526 outfile << "grp"<< dom1 <<".Add([ " << (int)sei+1 << " ])" << endl;
3527 if ( int dom2 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainOut())
3528 outfile << "grp"<< dom2 <<".Add([ " << (int)sei+1 << " ])" << endl;
3532 for (ElementIndex ei = 0; ei < ngMesh->GetNE(); ei++)
3534 Element el = (*ngMesh)[ei];
3535 outfile << "mesh.AddVolume([ ";
3536 for (int j = 0; j < el.GetNP(); j++)
3537 outfile << el[j] << ( j+1 < el.GetNP() ? ", " : " ])");
3541 for (int i = 1; i <= ngMesh->GetNSeg(); i++)
3543 const Segment & seg = ngMesh->LineSegment (i);
3544 outfile << "mesh.AddEdge([ "
3546 << seg[1] << " ])" << endl;
3548 cout << "Write " << pyFile << endl;
3551 //================================================================================
3553 * \brief Constructor of NETGENPlugin_ngMeshInfo
3555 //================================================================================
3557 NETGENPlugin_ngMeshInfo::NETGENPlugin_ngMeshInfo( netgen::Mesh* ngMesh):
3562 _nbNodes = ngMesh->GetNP();
3563 _nbSegments = ngMesh->GetNSeg();
3564 _nbFaces = ngMesh->GetNSE();
3565 _nbVolumes = ngMesh->GetNE();
3569 _nbNodes = _nbSegments = _nbFaces = _nbVolumes = 0;
3573 //================================================================================
3575 * \brief Copy LocalH member from one netgen mesh to another
3577 //================================================================================
3579 void NETGENPlugin_ngMeshInfo::transferLocalH( netgen::Mesh* fromMesh,
3580 netgen::Mesh* toMesh )
3582 if ( !fromMesh->LocalHFunctionGenerated() ) return;
3583 if ( !toMesh->LocalHFunctionGenerated() )
3585 toMesh->CalcLocalH(netgen::mparam.grading);
3587 toMesh->CalcLocalH();
3590 const size_t size = sizeof( netgen::LocalH );
3591 _copyOfLocalH = new char[ size ];
3592 memcpy( (void*)_copyOfLocalH, (void*)&toMesh->LocalHFunction(), size );
3593 memcpy( (void*)&toMesh->LocalHFunction(), (void*)&fromMesh->LocalHFunction(), size );
3596 //================================================================================
3598 * \brief Restore LocalH member of a netgen mesh
3600 //================================================================================
3602 void NETGENPlugin_ngMeshInfo::restoreLocalH( netgen::Mesh* toMesh )
3604 if ( _copyOfLocalH )
3606 const size_t size = sizeof( netgen::LocalH );
3607 memcpy( (void*)&toMesh->LocalHFunction(), (void*)_copyOfLocalH, size );
3608 delete [] _copyOfLocalH;
3613 //================================================================================
3615 * \brief Find "internal" sub-shapes
3617 //================================================================================
3619 NETGENPlugin_Internals::NETGENPlugin_Internals( SMESH_Mesh& mesh,
3620 const TopoDS_Shape& shape,
3622 : _mesh( mesh ), _is3D( is3D )
3624 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3626 TopExp_Explorer f,e;
3627 for ( f.Init( shape, TopAbs_FACE ); f.More(); f.Next() )
3629 int faceID = meshDS->ShapeToIndex( f.Current() );
3631 // find not computed internal edges
3633 for ( e.Init( f.Current().Oriented(TopAbs_FORWARD), TopAbs_EDGE ); e.More(); e.Next() )
3634 if ( e.Current().Orientation() == TopAbs_INTERNAL )
3636 SMESH_subMesh* eSM = mesh.GetSubMesh( e.Current() );
3637 if ( eSM->IsEmpty() )
3639 _e2face.insert( make_pair( eSM->GetId(), faceID ));
3640 for ( TopoDS_Iterator v(e.Current()); v.More(); v.Next() )
3641 _e2face.insert( make_pair( meshDS->ShapeToIndex( v.Value() ), faceID ));
3645 // find internal vertices in a face
3646 set<int> intVV; // issue 0020850 where same vertex is twice in a face
3647 for ( TopoDS_Iterator fSub( f.Current() ); fSub.More(); fSub.Next())
3648 if ( fSub.Value().ShapeType() == TopAbs_VERTEX )
3650 int vID = meshDS->ShapeToIndex( fSub.Value() );
3651 if ( intVV.insert( vID ).second )
3652 _f2v[ faceID ].push_back( vID );
3657 // find internal faces and their subshapes where nodes are to be doubled
3658 // to make a crack with non-sewed borders
3660 if ( f.Current().Orientation() == TopAbs_INTERNAL )
3662 _intShapes.insert( meshDS->ShapeToIndex( f.Current() ));
3665 list< TopoDS_Shape > edges;
3666 for ( e.Init( f.Current(), TopAbs_EDGE ); e.More(); e.Next())
3667 if ( SMESH_MesherHelper::NbAncestors( e.Current(), mesh, TopAbs_FACE ) > 1 )
3669 _intShapes.insert( meshDS->ShapeToIndex( e.Current() ));
3670 edges.push_back( e.Current() );
3671 // find border faces
3672 PShapeIteratorPtr fIt =
3673 SMESH_MesherHelper::GetAncestors( edges.back(),mesh,TopAbs_FACE );
3674 while ( const TopoDS_Shape* pFace = fIt->next() )
3675 if ( !pFace->IsSame( f.Current() ))
3676 _borderFaces.insert( meshDS->ShapeToIndex( *pFace ));
3679 // we consider vertex internal if it is shared by more than one internal edge
3680 list< TopoDS_Shape >::iterator edge = edges.begin();
3681 for ( ; edge != edges.end(); ++edge )
3682 for ( TopoDS_Iterator v( *edge ); v.More(); v.Next() )
3684 set<int> internalEdges;
3685 PShapeIteratorPtr eIt =
3686 SMESH_MesherHelper::GetAncestors( v.Value(),mesh,TopAbs_EDGE );
3687 while ( const TopoDS_Shape* pEdge = eIt->next() )
3689 int edgeID = meshDS->ShapeToIndex( *pEdge );
3690 if ( isInternalShape( edgeID ))
3691 internalEdges.insert( edgeID );
3693 if ( internalEdges.size() > 1 )
3694 _intShapes.insert( meshDS->ShapeToIndex( v.Value() ));
3698 } // loop on geom faces
3700 // find vertices internal in solids
3703 for ( TopExp_Explorer so(shape, TopAbs_SOLID); so.More(); so.Next())
3705 int soID = meshDS->ShapeToIndex( so.Current() );
3706 for ( TopoDS_Iterator soSub( so.Current() ); soSub.More(); soSub.Next())
3707 if ( soSub.Value().ShapeType() == TopAbs_VERTEX )
3708 _s2v[ soID ].push_back( meshDS->ShapeToIndex( soSub.Value() ));
3713 //================================================================================
3715 * \brief Find mesh faces on non-internal geom faces sharing internal edge
3716 * some nodes of which are to be doubled to make the second border of the "crack"
3718 //================================================================================
3720 void NETGENPlugin_Internals::findBorderElements( TIDSortedElemSet & borderElems )
3722 if ( _intShapes.empty() ) return;
3724 SMESH_Mesh& mesh = const_cast<SMESH_Mesh&>(_mesh);
3725 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3727 // loop on internal geom edges
3728 set<int>::const_iterator intShapeId = _intShapes.begin();
3729 for ( ; intShapeId != _intShapes.end(); ++intShapeId )
3731 const TopoDS_Shape& s = meshDS->IndexToShape( *intShapeId );
3732 if ( s.ShapeType() != TopAbs_EDGE ) continue;
3734 // get internal and non-internal geom faces sharing the internal edge <s>
3736 set<int>::iterator bordFace = _borderFaces.end();
3737 PShapeIteratorPtr faces = SMESH_MesherHelper::GetAncestors( s, _mesh, TopAbs_FACE );
3738 while ( const TopoDS_Shape* pFace = faces->next() )
3740 int faceID = meshDS->ShapeToIndex( *pFace );
3741 if ( isInternalShape( faceID ))
3744 bordFace = _borderFaces.insert( faceID ).first;
3746 if ( bordFace == _borderFaces.end() || !intFace ) continue;
3748 // get all links of mesh faces on internal geom face sharing nodes on edge <s>
3749 set< SMESH_OrientedLink > links; //!< links of faces on internal geom face
3750 list<const SMDS_MeshElement*> suspectFaces[2]; //!< mesh faces on border geom faces
3751 int nbSuspectFaces = 0;
3752 SMESHDS_SubMesh* intFaceSM = meshDS->MeshElements( intFace );
3753 if ( !intFaceSM || intFaceSM->NbElements() == 0 ) continue;
3754 SMESH_subMeshIteratorPtr smIt = mesh.GetSubMesh( s )->getDependsOnIterator(true,true);
3755 while ( smIt->more() )
3757 SMESHDS_SubMesh* sm = smIt->next()->GetSubMeshDS();
3758 if ( !sm ) continue;
3759 SMDS_NodeIteratorPtr nIt = sm->GetNodes();
3760 while ( nIt->more() )
3762 const SMDS_MeshNode* nOnEdge = nIt->next();
3763 SMDS_ElemIteratorPtr fIt = nOnEdge->GetInverseElementIterator(SMDSAbs_Face);
3764 while ( fIt->more() )
3766 const SMDS_MeshElement* f = fIt->next();
3767 const int nbNodes = f->NbCornerNodes();
3768 if ( intFaceSM->Contains( f ))
3770 for ( int i = 0; i < nbNodes; ++i )
3771 links.insert( SMESH_OrientedLink( f->GetNode(i), f->GetNode((i+1)%nbNodes)));
3776 for ( int i = 0; i < nbNodes; ++i )
3777 nbDblNodes += isInternalShape( f->GetNode(i)->getshapeId() );
3779 suspectFaces[ nbDblNodes < 2 ].push_back( f );
3785 // suspectFaces[0] having link with same orientation as mesh faces on
3786 // the internal geom face are <borderElems>. suspectFaces[1] have
3787 // only one node on edge <s>, we decide on them later (at the 2nd loop)
3788 // by links of <borderElems> found at the 1st and 2nd loops
3789 set< SMESH_OrientedLink > borderLinks;
3790 for ( int isPostponed = 0; isPostponed < 2; ++isPostponed )
3792 list<const SMDS_MeshElement*>::iterator fIt = suspectFaces[isPostponed].begin();
3793 for ( int nbF = 0; fIt != suspectFaces[isPostponed].end(); ++fIt, ++nbF )
3795 const SMDS_MeshElement* f = *fIt;
3796 bool isBorder = false, linkFound = false, borderLinkFound = false;
3797 list< SMESH_OrientedLink > faceLinks;
3798 int nbNodes = f->NbCornerNodes();
3799 for ( int i = 0; i < nbNodes; ++i )
3801 SMESH_OrientedLink link( f->GetNode(i), f->GetNode((i+1)%nbNodes));
3802 faceLinks.push_back( link );
3805 set< SMESH_OrientedLink >::iterator foundLink = links.find( link );
3806 if ( foundLink != links.end() )
3809 isBorder = ( foundLink->_reversed == link._reversed );
3810 if ( !isBorder && !isPostponed ) break;
3811 faceLinks.pop_back();
3813 else if ( isPostponed && !borderLinkFound )
3815 foundLink = borderLinks.find( link );
3816 if ( foundLink != borderLinks.end() )
3818 borderLinkFound = true;
3819 isBorder = ( foundLink->_reversed != link._reversed );
3826 borderElems.insert( f );
3827 borderLinks.insert( faceLinks.begin(), faceLinks.end() );
3829 else if ( !linkFound && !borderLinkFound )
3831 suspectFaces[1].push_back( f );
3832 if ( nbF > 2 * nbSuspectFaces )
3833 break; // dead loop protection
3840 //================================================================================
3842 * \brief put internal shapes in maps and fill in submeshes to precompute
3844 //================================================================================
3846 void NETGENPlugin_Internals::getInternalEdges( TopTools_IndexedMapOfShape& fmap,
3847 TopTools_IndexedMapOfShape& emap,
3848 TopTools_IndexedMapOfShape& vmap,
3849 list< SMESH_subMesh* > smToPrecompute[])
3851 if ( !hasInternalEdges() ) return;
3852 map<int,int>::const_iterator ev_face = _e2face.begin();
3853 for ( ; ev_face != _e2face.end(); ++ev_face )
3855 const TopoDS_Shape& ev = _mesh.GetMeshDS()->IndexToShape( ev_face->first );
3856 const TopoDS_Shape& face = _mesh.GetMeshDS()->IndexToShape( ev_face->second );
3858 ( ev.ShapeType() == TopAbs_EDGE ? emap : vmap ).Add( ev );
3860 //cout<<"INTERNAL EDGE or VERTEX "<<ev_face->first<<" on face "<<ev_face->second<<endl;
3862 smToPrecompute[ MeshDim_1D ].push_back( _mesh.GetSubMeshContaining( ev_face->first ));
3866 //================================================================================
3868 * \brief return shapes and submeshes to be meshed and already meshed boundary submeshes
3870 //================================================================================
3872 void NETGENPlugin_Internals::getInternalFaces( TopTools_IndexedMapOfShape& fmap,
3873 TopTools_IndexedMapOfShape& emap,
3874 list< SMESH_subMesh* >& intFaceSM,
3875 list< SMESH_subMesh* >& boundarySM)
3877 if ( !hasInternalFaces() ) return;
3879 // <fmap> and <emap> are for not yet meshed shapes
3880 // <intFaceSM> is for submeshes of faces
3881 // <boundarySM> is for meshed edges and vertices
3886 set<int> shapeIDs ( _intShapes );
3887 if ( !_borderFaces.empty() )
3888 shapeIDs.insert( _borderFaces.begin(), _borderFaces.end() );
3890 set<int>::const_iterator intS = shapeIDs.begin();
3891 for ( ; intS != shapeIDs.end(); ++intS )
3893 SMESH_subMesh* sm = _mesh.GetSubMeshContaining( *intS );
3895 if ( sm->GetSubShape().ShapeType() != TopAbs_FACE ) continue;
3897 intFaceSM.push_back( sm );
3899 // add submeshes of not computed internal faces
3900 if ( !sm->IsEmpty() ) continue;
3902 SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(true,true);
3903 while ( smIt->more() )
3906 const TopoDS_Shape& s = sm->GetSubShape();
3908 if ( sm->IsEmpty() )
3911 switch ( s.ShapeType() ) {
3912 case TopAbs_FACE: fmap.Add ( s ); break;
3913 case TopAbs_EDGE: emap.Add ( s ); break;
3919 if ( s.ShapeType() != TopAbs_FACE )
3920 boundarySM.push_back( sm );
3926 //================================================================================
3928 * \brief Return true if given shape is to be precomputed in order to be correctly
3929 * added to netgen mesh
3931 //================================================================================
3933 bool NETGENPlugin_Internals::isShapeToPrecompute(const TopoDS_Shape& s)
3935 int shapeID = _mesh.GetMeshDS()->ShapeToIndex( s );
3936 switch ( s.ShapeType() ) {
3937 case TopAbs_FACE : break; //return isInternalShape( shapeID ) || isBorderFace( shapeID );
3938 case TopAbs_EDGE : return isInternalEdge( shapeID );
3939 case TopAbs_VERTEX: break;
3945 //================================================================================
3947 * \brief Return SMESH
3949 //================================================================================
3951 SMESH_Mesh& NETGENPlugin_Internals::getMesh() const
3953 return const_cast<SMESH_Mesh&>( _mesh );
3956 //================================================================================
3958 * \brief Access to a counter of NETGENPlugin_NetgenLibWrapper instances
3960 //================================================================================
3962 int& NETGENPlugin_NetgenLibWrapper::instanceCounter()
3964 static int theCouner = 0;
3968 //================================================================================
3970 * \brief Initialize netgen library
3972 //================================================================================
3974 NETGENPlugin_NetgenLibWrapper::NETGENPlugin_NetgenLibWrapper()
3976 if ( instanceCounter() == 0 )
3979 ++instanceCounter();
3981 _isComputeOk = false;
3985 if ( !getenv( "KEEP_NETGEN_OUTPUT" ))
3987 // redirect all netgen output (mycout,myerr,cout) to _outputFileName
3988 _outputFileName = getOutputFileName();
3989 _ngcout = netgen::mycout;
3990 _ngcerr = netgen::myerr;
3991 netgen::mycout = new ofstream ( _outputFileName.c_str() );
3992 netgen::myerr = netgen::mycout;
3993 _coutBuffer = std::cout.rdbuf();
3995 cout << "NOTE: netgen output is redirected to file " << _outputFileName << endl;
3997 std::cout.rdbuf( netgen::mycout->rdbuf() );
4001 _ngMesh = Ng_NewMesh();
4004 //================================================================================
4006 * \brief Finish using netgen library
4008 //================================================================================
4010 NETGENPlugin_NetgenLibWrapper::~NETGENPlugin_NetgenLibWrapper()
4012 --instanceCounter();
4014 Ng_DeleteMesh( _ngMesh );
4018 std::cout.rdbuf( _coutBuffer );
4025 //================================================================================
4027 * \brief Set netgen mesh to delete at destruction
4029 //================================================================================
4031 void NETGENPlugin_NetgenLibWrapper::setMesh( Ng_Mesh* mesh )
4034 Ng_DeleteMesh( _ngMesh );
4038 //================================================================================
4040 * \brief Return a unique file name
4042 //================================================================================
4044 std::string NETGENPlugin_NetgenLibWrapper::getOutputFileName()
4046 std::string aTmpDir = SALOMEDS_Tool::GetTmpDir();
4048 TCollection_AsciiString aGenericName = (char*)aTmpDir.c_str();
4049 aGenericName += "NETGEN_";
4051 aGenericName += getpid();
4053 aGenericName += _getpid();
4055 aGenericName += "_";
4056 aGenericName += Abs((Standard_Integer)(long) aGenericName.ToCString());
4057 aGenericName += ".out";
4059 return aGenericName.ToCString();
4062 //================================================================================
4064 * \brief Remove "test.out" and "problemfaces" files in current directory
4066 //================================================================================
4068 void NETGENPlugin_NetgenLibWrapper::RemoveTmpFiles()
4070 bool rm = SMESH_File("test.out").remove() ;
4072 if ( rm && netgen::testout && instanceCounter() == 0 )
4074 delete netgen::testout;
4075 netgen::testout = 0;
4078 SMESH_File("problemfaces").remove();
4079 SMESH_File("occmesh.rep").remove();
4082 //================================================================================
4084 * \brief Remove file with netgen output
4086 //================================================================================
4088 void NETGENPlugin_NetgenLibWrapper::removeOutputFile()
4090 if ( !_outputFileName.empty() )
4094 delete netgen::mycout;
4095 netgen::mycout = _ngcout;
4096 netgen::myerr = _ngcerr;
4099 string tmpDir = SALOMEDS_Tool::GetDirFromPath ( _outputFileName );
4100 string aFileName = SALOMEDS_Tool::GetNameFromPath( _outputFileName ) + ".out";
4101 SALOMEDS::ListOfFileNames_var aFiles = new SALOMEDS::ListOfFileNames;
4103 aFiles[0] = aFileName.c_str();
4105 SALOMEDS_Tool::RemoveTemporaryFiles( tmpDir.c_str(), aFiles.in(), true );