1 // Copyright (C) 2007-2016 CEA/DEN, EDF R&D, OPEN CASCADE
3 // Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License, or (at your option) any later version.
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
23 // NETGENPlugin : C++ implementation
24 // File : NETGENPlugin_Mesher.cxx
25 // Author : Michael Sazonov (OCN)
28 //=============================================================================
30 #include "NETGENPlugin_Mesher.hxx"
31 #include "NETGENPlugin_Hypothesis_2D.hxx"
32 #include "NETGENPlugin_SimpleHypothesis_3D.hxx"
34 #include <SMDS_FaceOfNodes.hxx>
35 #include <SMDS_LinearEdge.hxx>
36 #include <SMDS_MeshElement.hxx>
37 #include <SMDS_MeshNode.hxx>
38 #include <SMESHDS_Mesh.hxx>
39 #include <SMESH_Block.hxx>
40 #include <SMESH_Comment.hxx>
41 #include <SMESH_ComputeError.hxx>
42 #include <SMESH_ControlPnt.hxx>
43 #include <SMESH_File.hxx>
44 #include <SMESH_Gen_i.hxx>
45 #include <SMESH_Mesh.hxx>
46 #include <SMESH_MesherHelper.hxx>
47 #include <SMESH_subMesh.hxx>
48 #include <StdMeshers_QuadToTriaAdaptor.hxx>
49 #include <StdMeshers_ViscousLayers2D.hxx>
51 #include <SALOMEDS_Tool.hxx>
53 #include <utilities.h>
55 #include <BRepBuilderAPI_Copy.hxx>
56 #include <BRep_Tool.hxx>
57 #include <Bnd_B3d.hxx>
58 #include <NCollection_Map.hxx>
59 #include <Standard_ErrorHandler.hxx>
60 #include <Standard_ProgramError.hxx>
61 #include <TColStd_MapOfInteger.hxx>
63 #include <TopExp_Explorer.hxx>
64 #include <TopTools_DataMapIteratorOfDataMapOfShapeInteger.hxx>
65 #include <TopTools_DataMapIteratorOfDataMapOfShapeShape.hxx>
66 #include <TopTools_DataMapOfShapeInteger.hxx>
67 #include <TopTools_DataMapOfShapeShape.hxx>
68 #include <TopTools_MapOfShape.hxx>
71 // Netgen include files
75 #include <occgeom.hpp>
76 #include <meshing.hpp>
77 //#include <ngexception.hpp>
80 extern int OCCGenerateMesh (OCCGeometry&, Mesh*&, MeshingParameters&, int, int);
82 extern int OCCGenerateMesh (OCCGeometry&, Mesh*&, int, int, char*);
84 //extern void OCCSetLocalMeshSize(OCCGeometry & geom, Mesh & mesh);
85 extern MeshingParameters mparam;
86 extern volatile multithreadt multithread;
87 extern bool merge_solids;
89 // values used for occgeo.facemeshstatus
90 enum EFaceMeshStatus { FACE_NOT_TREATED = 0,
102 using namespace nglib;
106 #define nodeVec_ACCESS(index) ((SMDS_MeshNode*) nodeVec.at((index)))
108 #define nodeVec_ACCESS(index) ((SMDS_MeshNode*) nodeVec[index])
111 #define NGPOINT_COORDS(p) p(0),p(1),p(2)
114 // dump elements added to ng mesh
115 //#define DUMP_SEGMENTS
116 //#define DUMP_TRIANGLES
117 //#define DUMP_TRIANGLES_SCRIPT "/tmp/trias.py" //!< debug AddIntVerticesInSolids()
120 TopTools_IndexedMapOfShape ShapesWithLocalSize;
121 std::map<int,double> VertexId2LocalSize;
122 std::map<int,double> EdgeId2LocalSize;
123 std::map<int,double> FaceId2LocalSize;
124 std::map<int,double> SolidId2LocalSize;
126 std::vector<SMESHUtils::ControlPnt> ControlPoints;
127 std::set<int> ShapesWithControlPoints; // <-- allows calling SetLocalSize() several times w/o recomputing ControlPoints
129 //=============================================================================
133 //=============================================================================
135 NETGENPlugin_Mesher::NETGENPlugin_Mesher (SMESH_Mesh* mesh,
136 const TopoDS_Shape& aShape,
142 _fineness(NETGENPlugin_Hypothesis::GetDefaultFineness()),
143 _isViscousLayers2D(false),
152 SetDefaultParameters();
153 ShapesWithLocalSize.Clear();
154 VertexId2LocalSize.clear();
155 EdgeId2LocalSize.clear();
156 FaceId2LocalSize.clear();
157 SolidId2LocalSize.clear();
158 ControlPoints.clear();
159 ShapesWithControlPoints.clear();
162 //================================================================================
166 //================================================================================
168 NETGENPlugin_Mesher::~NETGENPlugin_Mesher()
176 //================================================================================
178 * Set pointer to NETGENPlugin_Mesher* field of the holder, that will be
179 * nullified at destruction of this
181 //================================================================================
183 void NETGENPlugin_Mesher::SetSelfPointer( NETGENPlugin_Mesher ** ptr )
194 //================================================================================
196 * \brief Initialize global NETGEN parameters with default values
198 //================================================================================
200 void NETGENPlugin_Mesher::SetDefaultParameters()
202 netgen::MeshingParameters& mparams = netgen::mparam;
203 // maximal mesh edge size
204 mparams.maxh = 0;//NETGENPlugin_Hypothesis::GetDefaultMaxSize();
206 // minimal number of segments per edge
207 mparams.segmentsperedge = NETGENPlugin_Hypothesis::GetDefaultNbSegPerEdge();
208 // rate of growth of size between elements
209 mparams.grading = NETGENPlugin_Hypothesis::GetDefaultGrowthRate();
210 // safety factor for curvatures (elements per radius)
211 mparams.curvaturesafety = NETGENPlugin_Hypothesis::GetDefaultNbSegPerRadius();
212 // create elements of second order
213 mparams.secondorder = NETGENPlugin_Hypothesis::GetDefaultSecondOrder();
214 // quad-dominated surface meshing
218 mparams.quad = NETGENPlugin_Hypothesis_2D::GetDefaultQuadAllowed();
219 _fineness = NETGENPlugin_Hypothesis::GetDefaultFineness();
220 mparams.uselocalh = NETGENPlugin_Hypothesis::GetDefaultSurfaceCurvature();
221 netgen::merge_solids = NETGENPlugin_Hypothesis::GetDefaultFuseEdges();
224 //=============================================================================
228 //=============================================================================
230 void SetLocalSize(TopoDS_Shape GeomShape, double LocalSize)
232 if ( GeomShape.IsNull() ) return;
233 TopAbs_ShapeEnum GeomType = GeomShape.ShapeType();
234 if (GeomType == TopAbs_COMPOUND) {
235 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()) {
236 SetLocalSize(it.Value(), LocalSize);
241 if (! ShapesWithLocalSize.Contains(GeomShape))
242 key = ShapesWithLocalSize.Add(GeomShape);
244 key = ShapesWithLocalSize.FindIndex(GeomShape);
245 if (GeomType == TopAbs_VERTEX) {
246 VertexId2LocalSize[key] = LocalSize;
247 } else if (GeomType == TopAbs_EDGE) {
248 EdgeId2LocalSize[key] = LocalSize;
249 } else if (GeomType == TopAbs_FACE) {
250 FaceId2LocalSize[key] = LocalSize;
251 } else if (GeomType == TopAbs_SOLID) {
252 SolidId2LocalSize[key] = LocalSize;
256 //=============================================================================
258 * Pass parameters to NETGEN
260 //=============================================================================
261 void NETGENPlugin_Mesher::SetParameters(const NETGENPlugin_Hypothesis* hyp)
265 netgen::MeshingParameters& mparams = netgen::mparam;
266 // Initialize global NETGEN parameters:
267 // maximal mesh segment size
268 mparams.maxh = hyp->GetMaxSize();
269 // maximal mesh element linear size
270 mparams.minh = hyp->GetMinSize();
271 // minimal number of segments per edge
272 mparams.segmentsperedge = hyp->GetNbSegPerEdge();
273 // rate of growth of size between elements
274 mparams.grading = hyp->GetGrowthRate();
275 // safety factor for curvatures (elements per radius)
276 mparams.curvaturesafety = hyp->GetNbSegPerRadius();
277 // create elements of second order
278 mparams.secondorder = hyp->GetSecondOrder() ? 1 : 0;
279 // quad-dominated surface meshing
280 mparams.quad = hyp->GetQuadAllowed() ? 1 : 0;
281 _optimize = hyp->GetOptimize();
282 _fineness = hyp->GetFineness();
283 mparams.uselocalh = hyp->GetSurfaceCurvature();
284 netgen::merge_solids = hyp->GetFuseEdges();
287 SMESH_Gen_i* smeshGen_i = SMESH_Gen_i::GetSMESHGen();
288 CORBA::Object_var anObject = smeshGen_i->GetNS()->Resolve("/myStudyManager");
289 SALOMEDS::StudyManager_var aStudyMgr = SALOMEDS::StudyManager::_narrow(anObject);
290 SALOMEDS::Study_var myStudy = aStudyMgr->GetStudyByID(hyp->GetStudyId());
291 if ( !myStudy->_is_nil() )
293 const NETGENPlugin_Hypothesis::TLocalSize localSizes = hyp->GetLocalSizesAndEntries();
294 NETGENPlugin_Hypothesis::TLocalSize::const_iterator it = localSizes.begin();
295 for ( ; it != localSizes.end() ; it++)
297 std::string entry = (*it).first;
298 double val = (*it).second;
300 GEOM::GEOM_Object_var aGeomObj;
301 SALOMEDS::SObject_var aSObj = myStudy->FindObjectID( entry.c_str() );
302 if ( !aSObj->_is_nil() ) {
303 CORBA::Object_var obj = aSObj->GetObject();
304 aGeomObj = GEOM::GEOM_Object::_narrow(obj);
307 TopoDS_Shape S = smeshGen_i->GeomObjectToShape( aGeomObj.in() );
308 ::SetLocalSize(S, val);
314 //=============================================================================
316 * Pass simple parameters to NETGEN
318 //=============================================================================
320 void NETGENPlugin_Mesher::SetParameters(const NETGENPlugin_SimpleHypothesis_2D* hyp)
324 SetDefaultParameters();
327 //=============================================================================
329 * Link - a pair of integer numbers
331 //=============================================================================
335 Link(int _n1, int _n2) : n1(_n1), n2(_n2) {}
336 Link() : n1(0), n2(0) {}
337 bool Contains( int n ) const { return n == n1 || n == n2; }
338 bool IsConnected( const Link& other ) const
340 return (( Contains( other.n1 ) || Contains( other.n2 )) && ( this != &other ));
344 int HashCode(const Link& aLink, int aLimit)
346 return HashCode(aLink.n1 + aLink.n2, aLimit);
349 Standard_Boolean IsEqual(const Link& aLink1, const Link& aLink2)
351 return (( aLink1.n1 == aLink2.n1 && aLink1.n2 == aLink2.n2 ) ||
352 ( aLink1.n1 == aLink2.n2 && aLink1.n2 == aLink2.n1 ));
357 //================================================================================
359 * \brief return id of netgen point corresponding to SMDS node
361 //================================================================================
362 typedef map< const SMDS_MeshNode*, int > TNode2IdMap;
364 int ngNodeId( const SMDS_MeshNode* node,
365 netgen::Mesh& ngMesh,
366 TNode2IdMap& nodeNgIdMap)
368 int newNgId = ngMesh.GetNP() + 1;
370 TNode2IdMap::iterator node_id = nodeNgIdMap.insert( make_pair( node, newNgId )).first;
372 if ( node_id->second == newNgId)
374 #if defined(DUMP_SEGMENTS) || defined(DUMP_TRIANGLES)
375 cout << "Ng " << newNgId << " - " << node;
377 netgen::MeshPoint p( netgen::Point<3> (node->X(), node->Y(), node->Z()) );
378 ngMesh.AddPoint( p );
380 return node_id->second;
383 //================================================================================
385 * \brief Return computed EDGEs connected to the given one
387 //================================================================================
389 list< TopoDS_Edge > getConnectedEdges( const TopoDS_Edge& edge,
390 const TopoDS_Face& face,
391 const set< SMESH_subMesh* > & computedSM,
392 const SMESH_MesherHelper& helper,
393 map< SMESH_subMesh*, set< int > >& addedEdgeSM2Faces)
396 list< TopoDS_Edge > edges;
397 list< int > nbEdgesInWire;
398 /*int nbWires =*/ SMESH_Block::GetOrderedEdges( face, edges, nbEdgesInWire);
400 // find <edge> within <edges>
401 list< TopoDS_Edge >::iterator eItFwd = edges.begin();
402 for ( ; eItFwd != edges.end(); ++eItFwd )
403 if ( edge.IsSame( *eItFwd ))
405 if ( eItFwd == edges.end()) return list< TopoDS_Edge>();
407 if ( eItFwd->Orientation() >= TopAbs_INTERNAL )
409 // connected INTERNAL edges returned from GetOrderedEdges() are wrongly oriented
410 // so treat each INTERNAL edge separately
411 TopoDS_Edge e = *eItFwd;
413 edges.push_back( e );
417 // get all computed EDGEs connected to <edge>
419 list< TopoDS_Edge >::iterator eItBack = eItFwd, ePrev;
420 TopoDS_Vertex vCommon;
421 TopTools_MapOfShape eAdded; // map used not to add a seam edge twice to <edges>
424 // put edges before <edge> to <edges> back
425 while ( edges.begin() != eItFwd )
426 edges.splice( edges.end(), edges, edges.begin() );
430 while ( ++eItFwd != edges.end() )
432 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( *eItFwd );
434 bool connected = TopExp::CommonVertex( *ePrev, *eItFwd, vCommon );
435 bool computed = sm->IsMeshComputed();
436 bool added = addedEdgeSM2Faces[ sm ].count( helper.GetSubShapeID() );
437 bool doubled = !eAdded.Add( *eItFwd );
438 bool orientOK = (( ePrev ->Orientation() < TopAbs_INTERNAL ) ==
439 ( eItFwd->Orientation() < TopAbs_INTERNAL ) );
440 if ( !connected || !computed || !orientOK || added || doubled )
442 // stop advancement; move edges from tail to head
443 while ( edges.back() != *ePrev )
444 edges.splice( edges.begin(), edges, --edges.end() );
450 while ( eItBack != edges.begin() )
454 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( *eItBack );
456 bool connected = TopExp::CommonVertex( *ePrev, *eItBack, vCommon );
457 bool computed = sm->IsMeshComputed();
458 bool added = addedEdgeSM2Faces[ sm ].count( helper.GetSubShapeID() );
459 bool doubled = !eAdded.Add( *eItBack );
460 bool orientOK = (( ePrev ->Orientation() < TopAbs_INTERNAL ) ==
461 ( eItBack->Orientation() < TopAbs_INTERNAL ) );
462 if ( !connected || !computed || !orientOK || added || doubled)
465 edges.erase( edges.begin(), ePrev );
469 if ( edges.front() != edges.back() )
471 // assure that the 1st vertex is meshed
472 TopoDS_Edge eLast = edges.back();
473 while ( !SMESH_Algo::VertexNode( SMESH_MesherHelper::IthVertex( 0, edges.front()), helper.GetMeshDS())
475 edges.front() != eLast )
476 edges.splice( edges.end(), edges, edges.begin() );
481 //================================================================================
483 * \brief Make triangulation of a shape precise enough
485 //================================================================================
487 void updateTriangulation( const TopoDS_Shape& shape )
489 // static set< Poly_Triangulation* > updated;
491 // TopLoc_Location loc;
492 // TopExp_Explorer fExp( shape, TopAbs_FACE );
493 // for ( ; fExp.More(); fExp.Next() )
495 // Handle(Poly_Triangulation) triangulation =
496 // BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
497 // if ( triangulation.IsNull() ||
498 // updated.insert( triangulation.operator->() ).second )
500 // BRepTools::Clean (shape);
503 BRepMesh_IncrementalMesh e(shape, 0.01, true);
505 catch (Standard_Failure)
508 // updated.erase( triangulation.operator->() );
509 // triangulation = BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
510 // updated.insert( triangulation.operator->() );
514 //================================================================================
516 * \brief Returns a medium node either existing in SMESH of created by NETGEN
517 * \param [in] corner1 - corner node 1
518 * \param [in] corner2 - corner node 2
519 * \param [in] defaultMedium - the node created by NETGEN
520 * \param [in] helper - holder of medium nodes existing in SMESH
521 * \return const SMDS_MeshNode* - the result node
523 //================================================================================
525 const SMDS_MeshNode* mediumNode( const SMDS_MeshNode* corner1,
526 const SMDS_MeshNode* corner2,
527 const SMDS_MeshNode* defaultMedium,
528 const SMESH_MesherHelper* helper)
532 TLinkNodeMap::const_iterator l2n =
533 helper->GetTLinkNodeMap().find( SMESH_TLink( corner1, corner2 ));
534 if ( l2n != helper->GetTLinkNodeMap().end() )
535 defaultMedium = l2n->second;
537 return defaultMedium;
540 //================================================================================
542 * \brief Assure that mesh on given shapes is quadratic
544 //================================================================================
546 // void makeQuadratic( const TopTools_IndexedMapOfShape& shapes,
547 // SMESH_Mesh* mesh )
549 // for ( int i = 1; i <= shapes.Extent(); ++i )
551 // SMESHDS_SubMesh* smDS = mesh->GetMeshDS()->MeshElements( shapes(i) );
552 // if ( !smDS ) continue;
553 // SMDS_ElemIteratorPtr elemIt = smDS->GetElements();
554 // if ( !elemIt->more() ) continue;
555 // const SMDS_MeshElement* e = elemIt->next();
556 // if ( !e || e->IsQuadratic() )
559 // TIDSortedElemSet elems;
560 // elems.insert( e );
561 // while ( elemIt->more() )
562 // elems.insert( elems.end(), elemIt->next() );
564 // SMESH_MeshEditor( mesh ).ConvertToQuadratic( /*3d=*/false, elems, /*biQuad=*/false );
568 //================================================================================
570 * \brief Restrict size of elements on the given edge
572 //================================================================================
574 void setLocalSize(const TopoDS_Edge& edge,
578 if ( size <= std::numeric_limits<double>::min() )
580 Standard_Real u1, u2;
581 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, u1, u2);
582 if ( curve.IsNull() )
584 TopoDS_Iterator vIt( edge );
585 if ( !vIt.More() ) return;
586 gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex( vIt.Value() ));
587 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
591 const int nb = (int)( 1.5 * SMESH_Algo::EdgeLength( edge ) / size );
592 Standard_Real delta = (u2-u1)/nb;
593 for(int i=0; i<nb; i++)
595 Standard_Real u = u1 + delta*i;
596 gp_Pnt p = curve->Value(u);
597 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
598 netgen::Point3d pi(p.X(), p.Y(), p.Z());
599 double resultSize = mesh.GetH(pi);
600 if ( resultSize - size > 0.1*size )
601 // netgen does restriction iff oldH/newH > 1.2 (localh.cpp:136)
602 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), resultSize/1.201 );
608 //================================================================================
610 * \brief Set local size on shapes defined by SetParameters()
612 //================================================================================
614 void NETGENPlugin_Mesher::SetLocalSize( netgen::OCCGeometry& occgeo,
615 netgen::Mesh& ngMesh )
617 for(std::map<int,double>::const_iterator it=EdgeId2LocalSize.begin(); it!=EdgeId2LocalSize.end(); it++)
619 int key = (*it).first;
620 double hi = (*it).second;
621 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
622 setLocalSize( TopoDS::Edge(shape), hi, ngMesh );
624 for(std::map<int,double>::const_iterator it=VertexId2LocalSize.begin(); it!=VertexId2LocalSize.end(); it++)
626 int key = (*it).first;
627 double hi = (*it).second;
628 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
629 gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex(shape) );
630 NETGENPlugin_Mesher::RestrictLocalSize( ngMesh, p.XYZ(), hi );
632 for(map<int,double>::const_iterator it=FaceId2LocalSize.begin(); it!=FaceId2LocalSize.end(); it++)
634 int key = (*it).first;
635 double val = (*it).second;
636 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
637 int faceNgID = occgeo.fmap.FindIndex(shape);
640 occgeo.SetFaceMaxH(faceNgID, val);
641 for ( TopExp_Explorer edgeExp( shape, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
642 setLocalSize( TopoDS::Edge( edgeExp.Current() ), val, ngMesh );
644 else if ( !ShapesWithControlPoints.count( key ))
646 SMESHUtils::createPointsSampleFromFace( TopoDS::Face( shape ), val, ControlPoints );
647 ShapesWithControlPoints.insert( key );
650 for(map<int,double>::const_iterator it=SolidId2LocalSize.begin(); it!=SolidId2LocalSize.end(); it++)
652 int key = (*it).first;
653 double val = (*it).second;
654 if ( !ShapesWithControlPoints.count( key ))
656 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
657 SMESHUtils::createPointsSampleFromSolid( TopoDS::Solid( shape ), val, ControlPoints );
658 ShapesWithControlPoints.insert( key );
662 if ( !ControlPoints.empty() )
664 for ( size_t i = 1; i < ControlPoints.size(); ++i )
665 NETGENPlugin_Mesher::RestrictLocalSize( ngMesh, ControlPoints[i].XYZ(), ControlPoints[i].Size() );
669 //================================================================================
671 * \brief Initialize netgen::OCCGeometry with OCCT shape
673 //================================================================================
675 void NETGENPlugin_Mesher::PrepareOCCgeometry(netgen::OCCGeometry& occgeo,
676 const TopoDS_Shape& shape,
678 list< SMESH_subMesh* > * meshedSM,
679 NETGENPlugin_Internals* intern)
681 updateTriangulation( shape );
684 BRepBndLib::Add (shape, bb);
685 double x1,y1,z1,x2,y2,z2;
686 bb.Get (x1,y1,z1,x2,y2,z2);
687 netgen::Point<3> p1 = netgen::Point<3> (x1,y1,z1);
688 netgen::Point<3> p2 = netgen::Point<3> (x2,y2,z2);
689 occgeo.boundingbox = netgen::Box<3> (p1,p2);
691 occgeo.shape = shape;
694 // fill maps of shapes of occgeo with not yet meshed subshapes
696 // get root submeshes
697 list< SMESH_subMesh* > rootSM;
698 const int shapeID = mesh.GetMeshDS()->ShapeToIndex( shape );
699 if ( shapeID > 0 ) { // SMESH_subMesh with ID 0 may exist, don't use it!
700 rootSM.push_back( mesh.GetSubMesh( shape ));
703 for ( TopoDS_Iterator it( shape ); it.More(); it.Next() )
704 rootSM.push_back( mesh.GetSubMesh( it.Value() ));
709 // add subshapes of empty submeshes
710 list< SMESH_subMesh* >::iterator rootIt = rootSM.begin(), rootEnd = rootSM.end();
711 for ( ; rootIt != rootEnd; ++rootIt ) {
712 SMESH_subMesh * root = *rootIt;
713 SMESH_subMeshIteratorPtr smIt = root->getDependsOnIterator(/*includeSelf=*/true,
714 /*complexShapeFirst=*/true);
715 // to find a right orientation of subshapes (PAL20462)
716 TopTools_IndexedMapOfShape subShapes;
717 TopExp::MapShapes(root->GetSubShape(), subShapes);
718 while ( smIt->more() )
720 SMESH_subMesh* sm = smIt->next();
721 TopoDS_Shape shape = sm->GetSubShape();
722 totNbFaces += ( shape.ShapeType() == TopAbs_FACE );
723 if ( intern && intern->isShapeToPrecompute( shape ))
725 if ( !meshedSM || sm->IsEmpty() )
727 if ( shape.ShapeType() != TopAbs_VERTEX )
728 shape = subShapes( subShapes.FindIndex( shape ));// shape -> index -> oriented shape
729 if ( shape.Orientation() >= TopAbs_INTERNAL )
730 shape.Orientation( TopAbs_FORWARD ); // isuue 0020676
731 switch ( shape.ShapeType() ) {
732 case TopAbs_FACE : occgeo.fmap.Add( shape ); break;
733 case TopAbs_EDGE : occgeo.emap.Add( shape ); break;
734 case TopAbs_VERTEX: occgeo.vmap.Add( shape ); break;
735 case TopAbs_SOLID :occgeo.somap.Add( shape ); break;
739 // collect submeshes of meshed shapes
742 const int dim = SMESH_Gen::GetShapeDim( shape );
743 meshedSM[ dim ].push_back( sm );
747 occgeo.facemeshstatus.SetSize (totNbFaces);
748 occgeo.facemeshstatus = 0;
749 occgeo.face_maxh_modified.SetSize(totNbFaces);
750 occgeo.face_maxh_modified = 0;
751 occgeo.face_maxh.SetSize(totNbFaces);
752 occgeo.face_maxh = netgen::mparam.maxh;
755 //================================================================================
757 * \brief Return a default min size value suitable for the given geometry.
759 //================================================================================
761 double NETGENPlugin_Mesher::GetDefaultMinSize(const TopoDS_Shape& geom,
762 const double maxSize)
764 updateTriangulation( geom );
768 const int* pi[4] = { &i1, &i2, &i3, &i1 };
771 TopExp_Explorer fExp( geom, TopAbs_FACE );
772 for ( ; fExp.More(); fExp.Next() )
774 Handle(Poly_Triangulation) triangulation =
775 BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
776 if ( triangulation.IsNull() ) continue;
777 const double fTol = BRep_Tool::Tolerance( TopoDS::Face( fExp.Current() ));
778 const TColgp_Array1OfPnt& points = triangulation->Nodes();
779 const Poly_Array1OfTriangle& trias = triangulation->Triangles();
780 for ( int iT = trias.Lower(); iT <= trias.Upper(); ++iT )
782 trias(iT).Get( i1, i2, i3 );
783 for ( int j = 0; j < 3; ++j )
785 double dist2 = points(*pi[j]).SquareDistance( points( *pi[j+1] ));
786 if ( dist2 < minh && fTol*fTol < dist2 )
788 bb.Add( points(*pi[j]));
792 if ( minh > 0.25 * bb.SquareExtent() ) // simple geometry, rough triangulation
794 minh = 1e-3 * sqrt( bb.SquareExtent());
795 //cout << "BND BOX minh = " <<minh << endl;
799 minh = 3 * sqrt( minh ); // triangulation for visualization is rather fine
800 //cout << "TRIANGULATION minh = " <<minh << endl;
802 if ( minh > 0.5 * maxSize )
808 //================================================================================
810 * \brief Restrict size of elements at a given point
812 //================================================================================
814 void NETGENPlugin_Mesher::RestrictLocalSize(netgen::Mesh& ngMesh,
817 const bool overrideMinH)
819 if ( size <= std::numeric_limits<double>::min() )
821 if ( netgen::mparam.minh > size )
825 ngMesh.SetMinimalH( size );
826 netgen::mparam.minh = size;
830 size = netgen::mparam.minh;
833 netgen::Point3d pi(p.X(), p.Y(), p.Z());
834 ngMesh.RestrictLocalH( pi, size );
837 //================================================================================
839 * \brief fill ngMesh with nodes and elements of computed submeshes
841 //================================================================================
843 bool NETGENPlugin_Mesher::FillNgMesh(netgen::OCCGeometry& occgeom,
844 netgen::Mesh& ngMesh,
845 vector<const SMDS_MeshNode*>& nodeVec,
846 const list< SMESH_subMesh* > & meshedSM,
847 SMESH_MesherHelper* quadHelper,
848 SMESH_ProxyMesh::Ptr proxyMesh)
850 TNode2IdMap nodeNgIdMap;
851 for ( size_t i = 1; i < nodeVec.size(); ++i )
852 nodeNgIdMap.insert( make_pair( nodeVec[i], i ));
854 TopTools_MapOfShape visitedShapes;
855 map< SMESH_subMesh*, set< int > > visitedEdgeSM2Faces;
856 set< SMESH_subMesh* > computedSM( meshedSM.begin(), meshedSM.end() );
858 SMESH_MesherHelper helper (*_mesh);
860 int faceNgID = ngMesh.GetNFD();
862 list< SMESH_subMesh* >::const_iterator smIt, smEnd = meshedSM.end();
863 for ( smIt = meshedSM.begin(); smIt != smEnd; ++smIt )
865 SMESH_subMesh* sm = *smIt;
866 if ( !visitedShapes.Add( sm->GetSubShape() ))
869 const SMESHDS_SubMesh * smDS = sm->GetSubMeshDS();
870 if ( !smDS ) continue;
872 switch ( sm->GetSubShape().ShapeType() )
874 case TopAbs_EDGE: { // EDGE
875 // ----------------------
876 TopoDS_Edge geomEdge = TopoDS::Edge( sm->GetSubShape() );
877 if ( geomEdge.Orientation() >= TopAbs_INTERNAL )
878 geomEdge.Orientation( TopAbs_FORWARD ); // issue 0020676
880 // Add ng segments for each not meshed FACE the EDGE bounds
881 PShapeIteratorPtr fIt = helper.GetAncestors( geomEdge, *sm->GetFather(), TopAbs_FACE );
882 while ( const TopoDS_Shape * anc = fIt->next() )
884 faceNgID = occgeom.fmap.FindIndex( *anc );
886 continue; // meshed face
888 int faceSMDSId = helper.GetMeshDS()->ShapeToIndex( *anc );
889 if ( visitedEdgeSM2Faces[ sm ].count( faceSMDSId ))
890 continue; // already treated EDGE
892 TopoDS_Face face = TopoDS::Face( occgeom.fmap( faceNgID ));
893 if ( face.Orientation() >= TopAbs_INTERNAL )
894 face.Orientation( TopAbs_FORWARD ); // issue 0020676
896 // get all meshed EDGEs of the FACE connected to geomEdge (issue 0021140)
897 helper.SetSubShape( face );
898 list< TopoDS_Edge > edges = getConnectedEdges( geomEdge, face, computedSM, helper,
899 visitedEdgeSM2Faces );
901 continue; // wrong ancestor?
903 // find out orientation of <edges> within <face>
904 TopoDS_Edge eNotSeam = edges.front();
905 if ( helper.HasSeam() )
907 list< TopoDS_Edge >::iterator eIt = edges.begin();
908 while ( helper.IsRealSeam( *eIt )) ++eIt;
909 if ( eIt != edges.end() )
912 TopAbs_Orientation fOri = helper.GetSubShapeOri( face, eNotSeam );
913 bool isForwad = ( fOri == eNotSeam.Orientation() || fOri >= TopAbs_INTERNAL );
915 // get all nodes from connected <edges>
916 const bool isQuad = smDS->IsQuadratic();
917 StdMeshers_FaceSide fSide( face, edges, _mesh, isForwad, isQuad );
918 const vector<UVPtStruct>& points = fSide.GetUVPtStruct();
919 if ( points.empty() )
920 return false; // invalid node params?
921 int i, nbSeg = fSide.NbSegments();
923 // remember EDGEs of fSide to treat only once
924 for ( int iE = 0; iE < fSide.NbEdges(); ++iE )
925 visitedEdgeSM2Faces[ helper.GetMesh()->GetSubMesh( fSide.Edge(iE )) ].insert(faceSMDSId);
927 double otherSeamParam = 0;
932 int prevNgId = ngNodeId( points[0].node, ngMesh, nodeNgIdMap );
934 for ( i = 0; i < nbSeg; ++i )
936 const UVPtStruct& p1 = points[ i ];
937 const UVPtStruct& p2 = points[ i+1 ];
939 if ( p1.node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX ) //an EDGE begins
942 if ( helper.IsRealSeam( p1.node->getshapeId() ))
944 TopoDS_Edge e = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
945 isSeam = helper.IsRealSeam( e );
948 otherSeamParam = helper.GetOtherParam( helper.GetPeriodicIndex() & 1 ? p2.u : p2.v );
955 seg[1] = prevNgId = ngNodeId( p2.node, ngMesh, nodeNgIdMap );
956 // node param on curve
957 seg.epgeominfo[ 0 ].dist = p1.param;
958 seg.epgeominfo[ 1 ].dist = p2.param;
960 seg.epgeominfo[ 0 ].u = p1.u;
961 seg.epgeominfo[ 0 ].v = p1.v;
962 seg.epgeominfo[ 1 ].u = p2.u;
963 seg.epgeominfo[ 1 ].v = p2.v;
965 //geomEdge = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
966 //seg.epgeominfo[ 0 ].edgenr = seg.epgeominfo[ 1 ].edgenr = occgeom.emap.FindIndex( geomEdge );
968 //seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
969 seg.si = faceNgID; // = geom.fmap.FindIndex (face);
970 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
971 ngMesh.AddSegment (seg);
973 SMESH_TNodeXYZ np1( p1.node ), np2( p2.node );
974 RestrictLocalSize( ngMesh, 0.5*(np1+np2), (np1-np2).Modulus() );
977 cout << "Segment: " << seg.edgenr << " on SMESH face " << helper.GetMeshDS()->ShapeToIndex( face ) << endl
978 << "\tface index: " << seg.si << endl
979 << "\tp1: " << seg[0] << endl
980 << "\tp2: " << seg[1] << endl
981 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
982 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
983 //<< "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
984 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
985 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl;
986 //<< "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
990 if ( helper.GetPeriodicIndex() && 1 ) {
991 seg.epgeominfo[ 0 ].u = otherSeamParam;
992 seg.epgeominfo[ 1 ].u = otherSeamParam;
993 swap (seg.epgeominfo[0].v, seg.epgeominfo[1].v);
995 seg.epgeominfo[ 0 ].v = otherSeamParam;
996 seg.epgeominfo[ 1 ].v = otherSeamParam;
997 swap (seg.epgeominfo[0].u, seg.epgeominfo[1].u);
999 swap( seg[0], seg[1] );
1000 swap( seg.epgeominfo[0].dist, seg.epgeominfo[1].dist );
1001 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1002 ngMesh.AddSegment( seg );
1003 #ifdef DUMP_SEGMENTS
1004 cout << "Segment: " << seg.edgenr << endl
1005 << "\t is SEAM (reverse) of the previous. "
1006 << " Other " << (helper.GetPeriodicIndex() && 1 ? "U" : "V")
1007 << " = " << otherSeamParam << endl;
1010 else if ( fOri == TopAbs_INTERNAL )
1012 swap( seg[0], seg[1] );
1013 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1014 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1015 ngMesh.AddSegment( seg );
1016 #ifdef DUMP_SEGMENTS
1017 cout << "Segment: " << seg.edgenr << endl << "\t is REVERSE of the previous" << endl;
1021 } // loop on geomEdge ancestors
1023 if ( quadHelper ) // remember medium nodes of sub-meshes
1025 SMDS_ElemIteratorPtr edges = smDS->GetElements();
1026 while ( edges->more() )
1028 const SMDS_MeshElement* e = edges->next();
1029 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshEdge*>( e )))
1035 } // case TopAbs_EDGE
1037 case TopAbs_FACE: { // FACE
1038 // ----------------------
1039 const TopoDS_Face& geomFace = TopoDS::Face( sm->GetSubShape() );
1040 helper.SetSubShape( geomFace );
1041 bool isInternalFace = ( geomFace.Orientation() == TopAbs_INTERNAL );
1043 // Find solids the geomFace bounds
1044 int solidID1 = 0, solidID2 = 0;
1045 StdMeshers_QuadToTriaAdaptor* quadAdaptor =
1046 dynamic_cast<StdMeshers_QuadToTriaAdaptor*>( proxyMesh.get() );
1049 solidID1 = occgeom.somap.FindIndex( quadAdaptor->GetShape() );
1053 PShapeIteratorPtr solidIt = helper.GetAncestors( geomFace, *sm->GetFather(), TopAbs_SOLID);
1054 while ( const TopoDS_Shape * solid = solidIt->next() )
1056 int id = occgeom.somap.FindIndex ( *solid );
1057 if ( solidID1 && id != solidID1 ) solidID2 = id;
1061 // Add ng face descriptors of meshed faces
1063 ngMesh.AddFaceDescriptor( netgen::FaceDescriptor( faceNgID, solidID1, solidID2, 0 ));
1065 // if second oreder is required, even already meshed faces must be passed to NETGEN
1066 int fID = occgeom.fmap.Add( geomFace );
1067 if ( occgeom.facemeshstatus.Size() < fID ) occgeom.facemeshstatus.SetSize( fID );
1068 occgeom.facemeshstatus[ fID-1 ] = netgen::FACE_MESHED_OK;
1069 while ( fID < faceNgID ) // geomFace is already in occgeom.fmap, add a copy
1071 fID = occgeom.fmap.Add( BRepBuilderAPI_Copy( geomFace, /*copyGeom=*/false ));
1072 if ( occgeom.facemeshstatus.Size() < fID ) occgeom.facemeshstatus.SetSize( fID );
1073 occgeom.facemeshstatus[ fID-1 ] = netgen::FACE_MESHED_OK;
1075 // Problem with the second order in a quadrangular mesh remains.
1076 // 1) All quadrangles generated by NETGEN are moved to an inexistent face
1077 // by FillSMesh() (find "AddFaceDescriptor")
1078 // 2) Temporary triangles generated by StdMeshers_QuadToTriaAdaptor
1079 // are on faces where quadrangles were.
1080 // Due to these 2 points, wrong geom faces are used while conversion to quadratic
1081 // of the mentioned above quadrangles and triangles
1083 // Orient the face correctly in solidID1 (issue 0020206)
1084 bool reverse = false;
1086 TopoDS_Shape solid = occgeom.somap( solidID1 );
1087 TopAbs_Orientation faceOriInSolid = helper.GetSubShapeOri( solid, geomFace );
1088 if ( faceOriInSolid >= 0 )
1090 helper.IsReversedSubMesh( TopoDS::Face( geomFace.Oriented( faceOriInSolid )));
1093 // Add surface elements
1095 netgen::Element2d tri(3);
1096 tri.SetIndex( faceNgID );
1097 SMESH_TNodeXYZ xyz[3];
1099 #ifdef DUMP_TRIANGLES
1100 cout << "SMESH face " << helper.GetMeshDS()->ShapeToIndex( geomFace )
1101 << " internal="<<isInternalFace << endl;
1104 smDS = proxyMesh->GetSubMesh( geomFace );
1106 SMDS_ElemIteratorPtr faces = smDS->GetElements();
1107 while ( faces->more() )
1109 const SMDS_MeshElement* f = faces->next();
1110 if ( f->NbNodes() % 3 != 0 ) // not triangle
1112 PShapeIteratorPtr solidIt=helper.GetAncestors(geomFace,*sm->GetFather(),TopAbs_SOLID);
1113 if ( const TopoDS_Shape * solid = solidIt->next() )
1114 sm = _mesh->GetSubMesh( *solid );
1115 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1116 smError.reset( new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH,"Not triangle sub-mesh"));
1117 smError->myBadElements.push_back( f );
1121 for ( int i = 0; i < 3; ++i )
1123 const SMDS_MeshNode* node = f->GetNode( i ), * inFaceNode=0;
1126 // get node UV on face
1127 int shapeID = node->getshapeId();
1128 if ( helper.IsSeamShape( shapeID ))
1130 if ( helper.IsSeamShape( f->GetNodeWrap( i+1 )->getshapeId() ))
1131 inFaceNode = f->GetNodeWrap( i-1 );
1133 inFaceNode = f->GetNodeWrap( i+1 );
1135 gp_XY uv = helper.GetNodeUV( geomFace, node, inFaceNode );
1137 int ind = reverse ? 3-i : i+1;
1138 tri.GeomInfoPi(ind).u = uv.X();
1139 tri.GeomInfoPi(ind).v = uv.Y();
1140 tri.PNum (ind) = ngNodeId( node, ngMesh, nodeNgIdMap );
1143 // pass a triangle size to NG size-map
1144 double size = ( ( xyz[0] - xyz[1] ).Modulus() +
1145 ( xyz[1] - xyz[2] ).Modulus() +
1146 ( xyz[2] - xyz[0] ).Modulus() ) / 3;
1147 gp_XYZ gc = ( xyz[0] + xyz[1] + xyz[2] ) / 3;
1148 RestrictLocalSize( ngMesh, gc, size, /*overrideMinH=*/false );
1150 ngMesh.AddSurfaceElement (tri);
1151 #ifdef DUMP_TRIANGLES
1152 cout << tri << endl;
1155 if ( isInternalFace )
1157 swap( tri[1], tri[2] );
1158 ngMesh.AddSurfaceElement (tri);
1159 #ifdef DUMP_TRIANGLES
1160 cout << tri << endl;
1165 if ( quadHelper ) // remember medium nodes of sub-meshes
1167 SMDS_ElemIteratorPtr faces = smDS->GetElements();
1168 while ( faces->more() )
1170 const SMDS_MeshElement* f = faces->next();
1171 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshFace*>( f )))
1177 } // case TopAbs_FACE
1179 case TopAbs_VERTEX: { // VERTEX
1180 // --------------------------
1181 // issue 0021405. Add node only if a VERTEX is shared by a not meshed EDGE,
1182 // else netgen removes a free node and nodeVector becomes invalid
1183 PShapeIteratorPtr ansIt = helper.GetAncestors( sm->GetSubShape(),
1187 while ( const TopoDS_Shape* e = ansIt->next() )
1189 SMESH_subMesh* eSub = helper.GetMesh()->GetSubMesh( *e );
1190 if (( toAdd = ( eSub->IsEmpty() && !SMESH_Algo::isDegenerated( TopoDS::Edge( *e )))))
1195 SMDS_NodeIteratorPtr nodeIt = smDS->GetNodes();
1196 if ( nodeIt->more() )
1197 ngNodeId( nodeIt->next(), ngMesh, nodeNgIdMap );
1203 } // loop on submeshes
1206 nodeVec.resize( ngMesh.GetNP() + 1 );
1207 TNode2IdMap::iterator node_NgId, nodeNgIdEnd = nodeNgIdMap.end();
1208 for ( node_NgId = nodeNgIdMap.begin(); node_NgId != nodeNgIdEnd; ++node_NgId)
1209 nodeVec[ node_NgId->second ] = node_NgId->first;
1214 //================================================================================
1216 * \brief Duplicate mesh faces on internal geom faces
1218 //================================================================================
1220 void NETGENPlugin_Mesher::FixIntFaces(const netgen::OCCGeometry& occgeom,
1221 netgen::Mesh& ngMesh,
1222 NETGENPlugin_Internals& internalShapes)
1224 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1226 // find ng indices of internal faces
1228 for ( int ngFaceID = 1; ngFaceID <= occgeom.fmap.Extent(); ++ngFaceID )
1230 int smeshID = meshDS->ShapeToIndex( occgeom.fmap( ngFaceID ));
1231 if ( internalShapes.isInternalShape( smeshID ))
1232 ngFaceIds.insert( ngFaceID );
1234 if ( !ngFaceIds.empty() )
1237 int i, nbFaces = ngMesh.GetNSE();
1238 for ( i = 1; i <= nbFaces; ++i)
1240 netgen::Element2d elem = ngMesh.SurfaceElement(i);
1241 if ( ngFaceIds.count( elem.GetIndex() ))
1243 swap( elem[1], elem[2] );
1244 ngMesh.AddSurfaceElement (elem);
1250 //================================================================================
1252 * \brief Tries to heal the mesh on a FACE. The FACE is supposed to be partially
1253 * meshed due to NETGEN failure
1254 * \param [in] occgeom - geometry
1255 * \param [in,out] ngMesh - the mesh to fix
1256 * \param [inout] faceID - ID of the FACE to fix the mesh on
1257 * \return bool - is mesh is or becomes OK
1259 //================================================================================
1261 bool NETGENPlugin_Mesher::FixFaceMesh(const netgen::OCCGeometry& occgeom,
1262 netgen::Mesh& ngMesh,
1265 // we address a case where the FACE is almost fully meshed except small holes
1266 // of usually triangular shape at FACE boundary (IPAL52861)
1268 // The case appeared to be not simple: holes only look triangular but
1269 // indeed are a self intersecting polygon. A reason of the bug was in coincident
1270 // NG points on a seam edge. But the code below is very nice, leave it for
1275 if ( occgeom.fmap.Extent() < faceID )
1277 //const TopoDS_Face& face = TopoDS::Face( occgeom.fmap( faceID ));
1279 // find free links on the FACE
1280 NCollection_Map<Link> linkMap;
1281 for ( int iF = 1; iF <= ngMesh.GetNSE(); ++iF )
1283 const netgen::Element2d& elem = ngMesh.SurfaceElement(iF);
1284 if ( faceID != elem.GetIndex() )
1286 int n0 = elem[ elem.GetNP() - 1 ];
1287 for ( int i = 0; i < elem.GetNP(); ++i )
1290 Link link( n0, n1 );
1291 if ( !linkMap.Add( link ))
1292 linkMap.Remove( link );
1296 // add/remove boundary links
1297 for ( int iSeg = 1; iSeg <= ngMesh.GetNSeg(); ++iSeg )
1299 const netgen::Segment& seg = ngMesh.LineSegment( iSeg );
1300 if ( seg.si != faceID ) // !edgeIDs.Contains( seg.edgenr ))
1302 Link link( seg[1], seg[0] ); // reverse!!!
1303 if ( !linkMap.Add( link ))
1304 linkMap.Remove( link );
1306 if ( linkMap.IsEmpty() )
1308 if ( linkMap.Extent() < 3 )
1311 // make triangles of the links
1313 netgen::Element2d tri(3);
1314 tri.SetIndex ( faceID );
1316 NCollection_Map<Link>::Iterator linkIt( linkMap );
1317 Link link1 = linkIt.Value();
1318 // look for a link connected to link1
1319 NCollection_Map<Link>::Iterator linkIt2 = linkIt;
1320 for ( linkIt2.Next(); linkIt2.More(); linkIt2.Next() )
1322 const Link& link2 = linkIt2.Value();
1323 if ( link2.IsConnected( link1 ))
1325 // look for a link connected to both link1 and link2
1326 NCollection_Map<Link>::Iterator linkIt3 = linkIt2;
1327 for ( linkIt3.Next(); linkIt3.More(); linkIt3.Next() )
1329 const Link& link3 = linkIt3.Value();
1330 if ( link3.IsConnected( link1 ) &&
1331 link3.IsConnected( link2 ) )
1336 tri[2] = ( link2.Contains( link1.n1 ) ? link2.n1 : link3.n1 );
1337 if ( tri[0] == tri[2] || tri[1] == tri[2] )
1339 ngMesh.AddSurfaceElement( tri );
1341 // prepare for the next tria search
1342 if ( linkMap.Extent() == 3 )
1344 linkMap.Remove( link3 );
1345 linkMap.Remove( link2 );
1347 linkMap.Remove( link1 );
1348 link1 = linkIt.Value();
1361 //================================================================================
1362 // define gp_XY_Subtracted pointer to function calling gp_XY::Subtracted(gp_XY)
1363 gp_XY_FunPtr(Subtracted);
1364 //gp_XY_FunPtr(Added);
1366 //================================================================================
1368 * \brief Evaluate distance between two 2d points along the surface
1370 //================================================================================
1372 double evalDist( const gp_XY& uv1,
1374 const Handle(Geom_Surface)& surf,
1375 const int stopHandler=-1)
1377 if ( stopHandler > 0 ) // continue recursion
1379 gp_XY mid = SMESH_MesherHelper::GetMiddleUV( surf, uv1, uv2 );
1380 return evalDist( uv1,mid, surf, stopHandler-1 ) + evalDist( mid,uv2, surf, stopHandler-1 );
1382 double dist3D = surf->Value( uv1.X(), uv1.Y() ).Distance( surf->Value( uv2.X(), uv2.Y() ));
1383 if ( stopHandler == 0 ) // stop recursion
1386 // start recursion if necessary
1387 double dist2D = SMESH_MesherHelper::ApplyIn2D(surf, uv1, uv2, gp_XY_Subtracted, 0).Modulus();
1388 if ( fabs( dist3D - dist2D ) < dist2D * 1e-10 )
1389 return dist3D; // equal parametrization of a planar surface
1391 return evalDist( uv1, uv2, surf, 3 ); // start recursion
1394 //================================================================================
1396 * \brief Data of vertex internal in geom face
1398 //================================================================================
1402 gp_XY uv; //!< UV in face parametric space
1403 int ngId; //!< ng id of corrsponding node
1404 gp_XY uvClose; //!< UV of closest boundary node
1405 int ngIdClose; //!< ng id of closest boundary node
1408 //================================================================================
1410 * \brief Data of vertex internal in solid
1412 //================================================================================
1416 int ngId; //!< ng id of corresponding node
1417 int ngIdClose; //!< ng id of closest 2d mesh element
1418 int ngIdCloseN; //!< ng id of closest node of the closest 2d mesh element
1421 inline double dist2(const netgen::MeshPoint& p1, const netgen::MeshPoint& p2)
1423 return gp_Pnt( NGPOINT_COORDS(p1)).SquareDistance( gp_Pnt( NGPOINT_COORDS(p2)));
1427 //================================================================================
1429 * \brief Make netgen take internal vertices in faces into account by adding
1430 * segments including internal vertices
1432 * This function works in supposition that 1D mesh is already computed in ngMesh
1434 //================================================================================
1436 void NETGENPlugin_Mesher::AddIntVerticesInFaces(const netgen::OCCGeometry& occgeom,
1437 netgen::Mesh& ngMesh,
1438 vector<const SMDS_MeshNode*>& nodeVec,
1439 NETGENPlugin_Internals& internalShapes)
1441 if ((int) nodeVec.size() < ngMesh.GetNP() )
1442 nodeVec.resize( ngMesh.GetNP(), 0 );
1444 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1445 SMESH_MesherHelper helper( internalShapes.getMesh() );
1447 const map<int,list<int> >& face2Vert = internalShapes.getFacesWithVertices();
1448 map<int,list<int> >::const_iterator f2v = face2Vert.begin();
1449 for ( ; f2v != face2Vert.end(); ++f2v )
1451 const TopoDS_Face& face = TopoDS::Face( meshDS->IndexToShape( f2v->first ));
1452 if ( face.IsNull() ) continue;
1453 int faceNgID = occgeom.fmap.FindIndex (face);
1454 if ( faceNgID < 0 ) continue;
1456 TopLoc_Location loc;
1457 Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc);
1459 helper.SetSubShape( face );
1460 helper.SetElementsOnShape( true );
1462 // Get data of internal vertices and add them to ngMesh
1464 multimap< double, TIntVData > dist2VData; // sort vertices by distance from boundary nodes
1466 int i, nbSegInit = ngMesh.GetNSeg();
1468 // boundary characteristics
1469 double totSegLen2D = 0;
1472 const list<int>& iVertices = f2v->second;
1473 list<int>::const_iterator iv = iVertices.begin();
1474 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1477 // get node on vertex
1478 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1479 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1482 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1483 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1484 nV = SMESH_Algo::VertexNode( V, meshDS );
1485 if ( !nV ) continue;
1488 netgen::MeshPoint mp( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1489 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1490 vData.ngId = ngMesh.GetNP();
1491 nodeVec.push_back( nV );
1495 vData.uv = helper.GetNodeUV( face, nV, 0, &uvOK );
1496 if ( !uvOK ) helper.CheckNodeUV( face, nV, vData.uv, BRep_Tool::Tolerance(V),/*force=*/1);
1498 // loop on all segments of the face to find the node closest to vertex and to count
1499 // average segment 2d length
1500 double closeDist2 = numeric_limits<double>::max(), dist2;
1502 for (i = 1; i <= ngMesh.GetNSeg(); ++i)
1504 netgen::Segment & seg = ngMesh.LineSegment(i);
1505 if ( seg.si != faceNgID ) continue;
1507 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1509 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1510 if ( ngIdLast == seg[ iEnd ] ) continue;
1511 dist2 = helper.ApplyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1512 if ( dist2 < closeDist2 )
1513 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1514 ngIdLast = seg[ iEnd ];
1518 totSegLen2D += helper.ApplyIn2D(surf, uv[0], uv[1], gp_XY_Subtracted, false).Modulus();
1522 dist2VData.insert( make_pair( closeDist2, vData ));
1525 if ( totNbSeg == 0 ) break;
1526 double avgSegLen2d = totSegLen2D / totNbSeg;
1528 // Loop on vertices to add segments
1530 multimap< double, TIntVData >::iterator dist_vData = dist2VData.begin();
1531 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1533 double closeDist2 = dist_vData->first, dist2;
1534 TIntVData & vData = dist_vData->second;
1536 // try to find more close node among segments added for internal vertices
1537 for (i = nbSegInit+1; i <= ngMesh.GetNSeg(); ++i)
1539 netgen::Segment & seg = ngMesh.LineSegment(i);
1540 if ( seg.si != faceNgID ) continue;
1542 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1544 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1545 dist2 = helper.ApplyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1546 if ( dist2 < closeDist2 )
1547 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1550 // decide whether to use the closest node as the second end of segment or to
1551 // create a new point
1552 int segEnd1 = vData.ngId;
1553 int segEnd2 = vData.ngIdClose; // to use closest node
1554 gp_XY uvV = vData.uv, uvP = vData.uvClose;
1555 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1556 double nodeDist2D = sqrt( closeDist2 );
1557 double nodeDist3D = evalDist( vData.uv, vData.uvClose, surf );
1558 bool avgLenOK = ( avgSegLen2d < 0.75 * nodeDist2D );
1559 bool hintLenOK = ( segLenHint < 0.75 * nodeDist3D );
1560 //cout << "uvV " << uvV.X() <<","<<uvV.Y() << " ";
1561 if ( hintLenOK || avgLenOK )
1563 // create a point between the closest node and V
1566 double r = min( 0.5, ( hintLenOK ? segLenHint/nodeDist3D : avgSegLen2d/nodeDist2D ));
1567 // direction from V to closet node in 2D
1568 gp_Dir2d v2n( helper.ApplyIn2D(surf, uvP, uvV, gp_XY_Subtracted, false ));
1570 uvP = vData.uv + r * nodeDist2D * v2n.XY();
1571 gp_Pnt P = surf->Value( uvP.X(), uvP.Y() ).Transformed( loc );
1573 netgen::MeshPoint mp( netgen::Point<3> (P.X(), P.Y(), P.Z()));
1574 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1575 segEnd2 = ngMesh.GetNP();
1576 //cout << "Middle " << r << " uv " << uvP.X() << "," << uvP.Y() << "( " << ngMesh.Point(segEnd2).X()<<","<<ngMesh.Point(segEnd2).Y()<<","<<ngMesh.Point(segEnd2).Z()<<" )"<< endl;
1577 SMDS_MeshNode * nP = helper.AddNode(P.X(), P.Y(), P.Z());
1578 nodeVec.push_back( nP );
1580 //else cout << "at Node " << " uv " << uvP.X() << "," << uvP.Y() << endl;
1583 netgen::Segment seg;
1585 if ( segEnd1 > segEnd2 ) swap( segEnd1, segEnd2 ), swap( uvV, uvP );
1586 seg[0] = segEnd1; // ng node id
1587 seg[1] = segEnd2; // ng node id
1588 seg.edgenr = ngMesh.GetNSeg() + 1;// segment id
1591 seg.epgeominfo[ 0 ].dist = 0; // param on curve
1592 seg.epgeominfo[ 0 ].u = uvV.X();
1593 seg.epgeominfo[ 0 ].v = uvV.Y();
1594 seg.epgeominfo[ 1 ].dist = 1; // param on curve
1595 seg.epgeominfo[ 1 ].u = uvP.X();
1596 seg.epgeominfo[ 1 ].v = uvP.Y();
1598 // seg.epgeominfo[ 0 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1599 // seg.epgeominfo[ 1 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1601 ngMesh.AddSegment (seg);
1603 // add reverse segment
1604 swap( seg[0], seg[1] );
1605 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1606 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1607 ngMesh.AddSegment (seg);
1613 //================================================================================
1615 * \brief Make netgen take internal vertices in solids into account by adding
1616 * faces including internal vertices
1618 * This function works in supposition that 2D mesh is already computed in ngMesh
1620 //================================================================================
1622 void NETGENPlugin_Mesher::AddIntVerticesInSolids(const netgen::OCCGeometry& occgeom,
1623 netgen::Mesh& ngMesh,
1624 vector<const SMDS_MeshNode*>& nodeVec,
1625 NETGENPlugin_Internals& internalShapes)
1627 #ifdef DUMP_TRIANGLES_SCRIPT
1628 // create a python script making a mesh containing triangles added for internal vertices
1629 ofstream py(DUMP_TRIANGLES_SCRIPT);
1630 py << "import SMESH"<< endl
1631 << "from salome.smesh import smeshBuilder"<<endl
1632 << "smesh = smeshBuilder.New(salome.myStudy)"<<endl
1633 << "m = smesh.Mesh(name='triangles')" << endl;
1635 if ((int) nodeVec.size() < ngMesh.GetNP() )
1636 nodeVec.resize( ngMesh.GetNP(), 0 );
1638 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1639 SMESH_MesherHelper helper( internalShapes.getMesh() );
1641 const map<int,list<int> >& so2Vert = internalShapes.getSolidsWithVertices();
1642 map<int,list<int> >::const_iterator s2v = so2Vert.begin();
1643 for ( ; s2v != so2Vert.end(); ++s2v )
1645 const TopoDS_Shape& solid = meshDS->IndexToShape( s2v->first );
1646 if ( solid.IsNull() ) continue;
1647 int solidNgID = occgeom.somap.FindIndex (solid);
1648 if ( solidNgID < 0 && !occgeom.somap.IsEmpty() ) continue;
1650 helper.SetSubShape( solid );
1651 helper.SetElementsOnShape( true );
1653 // find ng indices of faces within the solid
1655 for (TopExp_Explorer fExp(solid, TopAbs_FACE); fExp.More(); fExp.Next() )
1656 ngFaceIds.insert( occgeom.fmap.FindIndex( fExp.Current() ));
1657 if ( ngFaceIds.size() == 1 && *ngFaceIds.begin() == 0 )
1658 ngFaceIds.insert( 1 );
1660 // Get data of internal vertices and add them to ngMesh
1662 multimap< double, TIntVSoData > dist2VData; // sort vertices by distance from ng faces
1664 int i, nbFaceInit = ngMesh.GetNSE();
1666 // boundary characteristics
1667 double totSegLen = 0;
1670 const list<int>& iVertices = s2v->second;
1671 list<int>::const_iterator iv = iVertices.begin();
1672 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1675 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1677 // get node on vertex
1678 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1681 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1682 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1683 nV = SMESH_Algo::VertexNode( V, meshDS );
1684 if ( !nV ) continue;
1687 netgen::MeshPoint mpV( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1688 ngMesh.AddPoint ( mpV, 1, netgen::FIXEDPOINT );
1689 vData.ngId = ngMesh.GetNP();
1690 nodeVec.push_back( nV );
1692 // loop on all 2d elements to find the one closest to vertex and to count
1693 // average segment length
1694 double closeDist2 = numeric_limits<double>::max(), avgDist2;
1695 for (i = 1; i <= ngMesh.GetNSE(); ++i)
1697 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1698 if ( !ngFaceIds.count( elem.GetIndex() )) continue;
1700 multimap< double, int> dist2nID; // sort nodes of element by distance from V
1701 for ( int j = 0; j < elem.GetNP(); ++j)
1703 netgen::MeshPoint mp = ngMesh.Point( elem[j] );
1704 double d2 = dist2( mpV, mp );
1705 dist2nID.insert( make_pair( d2, elem[j] ));
1706 avgDist2 += d2 / elem.GetNP();
1708 totNbSeg++, totSegLen+= sqrt( dist2( mp, ngMesh.Point( elem[(j+1)%elem.GetNP()])));
1710 double dist = dist2nID.begin()->first; //avgDist2;
1711 if ( dist < closeDist2 )
1712 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= dist;
1714 dist2VData.insert( make_pair( closeDist2, vData ));
1717 if ( totNbSeg == 0 ) break;
1718 double avgSegLen = totSegLen / totNbSeg;
1720 // Loop on vertices to add triangles
1722 multimap< double, TIntVSoData >::iterator dist_vData = dist2VData.begin();
1723 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1725 double closeDist2 = dist_vData->first;
1726 TIntVSoData & vData = dist_vData->second;
1728 const netgen::MeshPoint& mpV = ngMesh.Point( vData.ngId );
1730 // try to find more close face among ones added for internal vertices
1731 for (i = nbFaceInit+1; i <= ngMesh.GetNSE(); ++i)
1733 double avgDist2 = 0;
1734 multimap< double, int> dist2nID;
1735 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1736 for ( int j = 0; j < elem.GetNP(); ++j)
1738 double d = dist2( mpV, ngMesh.Point( elem[j] ));
1739 dist2nID.insert( make_pair( d, elem[j] ));
1740 avgDist2 += d / elem.GetNP();
1741 if ( avgDist2 < closeDist2 )
1742 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= avgDist2;
1745 // sort nodes of the closest face by angle with vector from V to the closest node
1746 const double tol = numeric_limits<double>::min();
1747 map< double, int > angle2ID;
1748 const netgen::Element2d& closeFace = ngMesh.SurfaceElement( vData.ngIdClose );
1749 netgen::MeshPoint mp[2];
1750 mp[0] = ngMesh.Point( vData.ngIdCloseN );
1751 gp_XYZ p1( NGPOINT_COORDS( mp[0] ));
1752 gp_XYZ pV( NGPOINT_COORDS( mpV ));
1753 gp_Vec v2p1( pV, p1 );
1754 double distN1 = v2p1.Magnitude();
1755 if ( distN1 <= tol ) continue;
1757 for ( int j = 0; j < closeFace.GetNP(); ++j)
1759 mp[1] = ngMesh.Point( closeFace[j] );
1760 gp_Vec v2p( pV, gp_Pnt( NGPOINT_COORDS( mp[1] )) );
1761 angle2ID.insert( make_pair( v2p1.Angle( v2p ), closeFace[j]));
1763 // get node with angle of 60 degrees or greater
1764 map< double, int >::iterator angle_id = angle2ID.lower_bound( 60. * M_PI / 180. );
1765 if ( angle_id == angle2ID.end() ) angle_id = --angle2ID.end();
1766 const double minAngle = 30. * M_PI / 180.;
1767 const double angle = angle_id->first;
1768 bool angleOK = ( angle > minAngle );
1770 // find points to create a triangle
1771 netgen::Element2d tri(3);
1773 tri[0] = vData.ngId;
1774 tri[1] = vData.ngIdCloseN; // to use the closest nodes
1775 tri[2] = angle_id->second; // to use the node with best angle
1777 // decide whether to use the closest node and the node with best angle or to create new ones
1778 for ( int isBestAngleN = 0; isBestAngleN < 2; ++isBestAngleN )
1780 bool createNew = !angleOK; //, distOK = true;
1782 int triInd = isBestAngleN ? 2 : 1;
1783 mp[isBestAngleN] = ngMesh.Point( tri[triInd] );
1788 double distN2 = sqrt( dist2( mpV, mp[isBestAngleN]));
1789 createNew = ( fabs( distN2 - distN1 ) > 0.25 * distN1 );
1791 else if ( angle < tol )
1793 v2p1.SetX( v2p1.X() + 1e-3 );
1799 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1800 bool avgLenOK = ( avgSegLen < 0.75 * distN1 );
1801 bool hintLenOK = ( segLenHint < 0.75 * distN1 );
1802 createNew = (createNew || avgLenOK || hintLenOK );
1803 // we create a new node not closer than 0.5 to the closest face
1804 // in order not to clash with other close face
1805 double r = min( 0.5, ( hintLenOK ? segLenHint : avgSegLen ) / distN1 );
1806 distFromV = r * distN1;
1810 // create a new point, between the node and the vertex if angleOK
1811 gp_XYZ p( NGPOINT_COORDS( mp[isBestAngleN] ));
1812 gp_Vec v2p( pV, p ); v2p.Normalize();
1813 if ( isBestAngleN && !angleOK )
1814 p = p1 + gp_Dir( v2p.XYZ() - v2p1.XYZ()).XYZ() * distN1 * 0.95;
1816 p = pV + v2p.XYZ() * distFromV;
1818 if ( !isBestAngleN ) p1 = p, distN1 = distFromV;
1820 mp[isBestAngleN].SetPoint( netgen::Point<3> (p.X(), p.Y(), p.Z()));
1821 ngMesh.AddPoint ( mp[isBestAngleN], 1, netgen::SURFACEPOINT );
1822 tri[triInd] = ngMesh.GetNP();
1823 nodeVec.push_back( helper.AddNode( p.X(), p.Y(), p.Z()) );
1826 ngMesh.AddSurfaceElement (tri);
1827 swap( tri[1], tri[2] );
1828 ngMesh.AddSurfaceElement (tri);
1830 #ifdef DUMP_TRIANGLES_SCRIPT
1831 py << "n1 = m.AddNode( "<< mpV(0)<<", "<< mpV(1)<<", "<< mpV(2)<<") "<< endl
1832 << "n2 = m.AddNode( "<< mp[0](0)<<", "<< mp[0](1)<<", "<< mp[0](2)<<") "<< endl
1833 << "n3 = m.AddNode( "<< mp[1](0)<<", "<< mp[1](1)<<", "<< mp[1](2)<<" )" << endl
1834 << "m.AddFace([n1,n2,n3])" << endl;
1836 } // loop on internal vertices of a solid
1838 } // loop on solids with internal vertices
1841 //================================================================================
1843 * \brief Fill netgen mesh with segments of a FACE
1844 * \param ngMesh - netgen mesh
1845 * \param geom - container of OCCT geometry to mesh
1846 * \param wires - data of nodes on FACE boundary
1847 * \param helper - mesher helper holding the FACE
1848 * \param nodeVec - vector of nodes in which node index == netgen ID
1849 * \retval SMESH_ComputeErrorPtr - error description
1851 //================================================================================
1853 SMESH_ComputeErrorPtr
1854 NETGENPlugin_Mesher::AddSegmentsToMesh(netgen::Mesh& ngMesh,
1855 netgen::OCCGeometry& geom,
1856 const TSideVector& wires,
1857 SMESH_MesherHelper& helper,
1858 vector< const SMDS_MeshNode* > & nodeVec,
1859 const bool overrideMinH)
1861 // ----------------------------
1862 // Check wires and count nodes
1863 // ----------------------------
1865 for ( size_t iW = 0; iW < wires.size(); ++iW )
1867 StdMeshers_FaceSidePtr wire = wires[ iW ];
1868 if ( wire->MissVertexNode() )
1870 // Commented for issue 0020960. It worked for the case, let's wait for case where it doesn't.
1871 // It seems that there is no reason for this limitation
1873 // (new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH, "Missing nodes on vertices"));
1875 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1876 if ((int) uvPtVec.size() != wire->NbPoints() )
1877 return SMESH_ComputeError::New(COMPERR_BAD_INPUT_MESH,
1878 SMESH_Comment("Unexpected nb of points on wire ") << iW
1879 << ": " << uvPtVec.size()<<" != "<<wire->NbPoints());
1880 nbNodes += wire->NbPoints();
1882 nodeVec.reserve( nodeVec.size() + nbNodes + 1 );
1883 if ( nodeVec.empty() )
1884 nodeVec.push_back( 0 );
1886 // -----------------
1888 // -----------------
1890 const bool wasNgMeshEmpty = ( ngMesh.GetNP() < 1 ); /* true => this method is called by
1891 NETGENPlugin_NETGEN_2D_ONLY */
1893 // map for nodes on vertices since they can be shared between wires
1894 // ( issue 0020676, face_int_box.brep) and nodes built by NETGEN
1895 map<const SMDS_MeshNode*, int > node2ngID;
1896 if ( !wasNgMeshEmpty ) // fill node2ngID with nodes built by NETGEN
1898 set< int > subIDs; // ids of sub-shapes of the FACE
1899 for ( size_t iW = 0; iW < wires.size(); ++iW )
1901 StdMeshers_FaceSidePtr wire = wires[ iW ];
1902 for ( int iE = 0, nbE = wire->NbEdges(); iE < nbE; ++iE )
1904 subIDs.insert( wire->EdgeID( iE ));
1905 subIDs.insert( helper.GetMeshDS()->ShapeToIndex( wire->FirstVertex( iE )));
1908 for ( size_t ngID = 1; ngID < nodeVec.size(); ++ngID )
1909 if ( subIDs.count( nodeVec[ngID]->getshapeId() ))
1910 node2ngID.insert( make_pair( nodeVec[ngID], ngID ));
1913 const int solidID = 0, faceID = geom.fmap.FindIndex( helper.GetSubShape() );
1914 if ( ngMesh.GetNFD() < 1 )
1915 ngMesh.AddFaceDescriptor( netgen::FaceDescriptor( faceID, solidID, solidID, 0 ));
1917 for ( size_t iW = 0; iW < wires.size(); ++iW )
1919 StdMeshers_FaceSidePtr wire = wires[ iW ];
1920 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1921 const int nbSegments = wire->NbPoints() - 1;
1923 // assure the 1st node to be in node2ngID, which is needed to correctly
1924 // "close chain of segments" (see below) in case if the 1st node is not
1925 // onVertex because it is on a Viscous layer
1926 node2ngID.insert( make_pair( uvPtVec[ 0 ].node, ngMesh.GetNP() + 1 ));
1928 // compute length of every segment
1929 vector<double> segLen( nbSegments );
1930 for ( int i = 0; i < nbSegments; ++i )
1931 segLen[i] = SMESH_TNodeXYZ( uvPtVec[ i ].node ).Distance( uvPtVec[ i+1 ].node );
1933 int edgeID = 1, posID = -2;
1934 bool isInternalWire = false;
1935 double vertexNormPar = 0;
1936 const int prevNbNGSeg = ngMesh.GetNSeg();
1937 for ( int i = 0; i < nbSegments; ++i ) // loop on segments
1939 // Add the first point of a segment
1941 const SMDS_MeshNode * n = uvPtVec[ i ].node;
1942 const int posShapeID = n->getshapeId();
1943 bool onVertex = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX );
1944 bool onEdge = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_EDGE );
1946 // skip nodes on degenerated edges
1947 if ( helper.IsDegenShape( posShapeID ) &&
1948 helper.IsDegenShape( uvPtVec[ i+1 ].node->getshapeId() ))
1951 int ngID1 = ngMesh.GetNP() + 1, ngID2 = ngID1+1;
1952 if ( onVertex || ( !wasNgMeshEmpty && onEdge ) || helper.IsRealSeam( posShapeID ))
1953 ngID1 = node2ngID.insert( make_pair( n, ngID1 )).first->second;
1954 if ( ngID1 > ngMesh.GetNP() )
1956 netgen::MeshPoint mp( netgen::Point<3> (n->X(), n->Y(), n->Z()) );
1957 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1958 nodeVec.push_back( n );
1960 else // n is in ngMesh already, and ngID2 in prev segment is wrong
1962 ngID2 = ngMesh.GetNP() + 1;
1963 if ( i > 0 ) // prev segment belongs to same wire
1965 netgen::Segment& prevSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
1972 netgen::Segment seg;
1974 seg[0] = ngID1; // ng node id
1975 seg[1] = ngID2; // ng node id
1976 seg.edgenr = ngMesh.GetNSeg() + 1; // ng segment id
1977 seg.si = faceID; // = geom.fmap.FindIndex (face);
1979 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1981 const UVPtStruct& pnt = uvPtVec[ i + iEnd ];
1983 seg.epgeominfo[ iEnd ].dist = pnt.param; // param on curve
1984 seg.epgeominfo[ iEnd ].u = pnt.u;
1985 seg.epgeominfo[ iEnd ].v = pnt.v;
1987 // find out edge id and node parameter on edge
1988 onVertex = ( pnt.normParam + 1e-10 > vertexNormPar );
1989 if ( onVertex || posShapeID != posID )
1992 double normParam = pnt.normParam;
1994 normParam = 0.5 * ( uvPtVec[ i ].normParam + uvPtVec[ i+1 ].normParam );
1995 int edgeIndexInWire = wire->EdgeIndex( normParam );
1996 vertexNormPar = wire->LastParameter( edgeIndexInWire );
1997 const TopoDS_Edge& edge = wire->Edge( edgeIndexInWire );
1998 edgeID = geom.emap.FindIndex( edge );
2000 isInternalWire = ( edge.Orientation() == TopAbs_INTERNAL );
2001 // if ( onVertex ) // param on curve is different on each of two edges
2002 // seg.epgeominfo[ iEnd ].dist = helper.GetNodeU( edge, pnt.node );
2004 seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
2007 ngMesh.AddSegment (seg);
2009 // restrict size of elements near the segment
2010 SMESH_TNodeXYZ np1( n ), np2( uvPtVec[ i+1 ].node );
2011 // get an average size of adjacent segments to avoid sharp change of
2012 // element size (regression on issue 0020452, note 0010898)
2013 int iPrev = SMESH_MesherHelper::WrapIndex( i-1, nbSegments );
2014 int iNext = SMESH_MesherHelper::WrapIndex( i+1, nbSegments );
2015 double sumH = segLen[ iPrev ] + segLen[ i ] + segLen[ iNext ];
2016 int nbSeg = ( int( segLen[ iPrev ] > sumH / 100.) +
2017 int( segLen[ i ] > sumH / 100.) +
2018 int( segLen[ iNext ] > sumH / 100.));
2020 RestrictLocalSize( ngMesh, 0.5*(np1+np2), sumH / nbSeg, overrideMinH );
2022 if ( isInternalWire )
2024 swap (seg[0], seg[1]);
2025 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
2026 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
2027 ngMesh.AddSegment (seg);
2029 } // loop on segments on a wire
2031 // close chain of segments
2032 if ( nbSegments > 0 )
2034 netgen::Segment& lastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() - int( isInternalWire ));
2035 const SMDS_MeshNode * lastNode = uvPtVec.back().node;
2036 lastSeg[1] = node2ngID.insert( make_pair( lastNode, lastSeg[1] )).first->second;
2037 if ( lastSeg[1] > ngMesh.GetNP() )
2039 netgen::MeshPoint mp( netgen::Point<3> (lastNode->X(), lastNode->Y(), lastNode->Z()) );
2040 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
2041 nodeVec.push_back( lastNode );
2043 if ( isInternalWire )
2045 netgen::Segment& realLastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
2046 realLastSeg[0] = lastSeg[1];
2050 #ifdef DUMP_SEGMENTS
2051 cout << "BEGIN WIRE " << iW << endl;
2052 for ( int i = prevNbNGSeg+1; i <= ngMesh.GetNSeg(); ++i )
2054 netgen::Segment& seg = ngMesh.LineSegment( i );
2056 netgen::Segment& prevSeg = ngMesh.LineSegment( i-1 );
2057 if ( seg[0] == prevSeg[1] && seg[1] == prevSeg[0] )
2059 cout << "Segment: " << seg.edgenr << endl << "\tis REVRESE of the previous one" << endl;
2063 cout << "Segment: " << seg.edgenr << endl
2064 << "\tp1: " << seg[0] << " n" << nodeVec[ seg[0]]->GetID() << endl
2065 << "\tp2: " << seg[1] << " n" << nodeVec[ seg[1]]->GetID() << endl
2066 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
2067 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
2068 << "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
2069 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
2070 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl
2071 << "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
2073 cout << "--END WIRE " << iW << endl;
2075 SMESH_Comment __not_unused_variable( prevNbNGSeg );
2078 } // loop on WIREs of a FACE
2080 // add a segment instead of an internal vertex
2081 if ( wasNgMeshEmpty )
2083 NETGENPlugin_Internals intShapes( *helper.GetMesh(), helper.GetSubShape(), /*is3D=*/false );
2084 AddIntVerticesInFaces( geom, ngMesh, nodeVec, intShapes );
2086 ngMesh.CalcSurfacesOfNode();
2091 //================================================================================
2093 * \brief Fill SMESH mesh according to contents of netgen mesh
2094 * \param occgeo - container of OCCT geometry to mesh
2095 * \param ngMesh - netgen mesh
2096 * \param initState - bn of entities in netgen mesh before computing
2097 * \param sMesh - SMESH mesh to fill in
2098 * \param nodeVec - vector of nodes in which node index == netgen ID
2099 * \param comment - returns problem description
2100 * \param quadHelper - holder of medium nodes of sub-meshes
2101 * \retval int - error
2103 //================================================================================
2105 int NETGENPlugin_Mesher::FillSMesh(const netgen::OCCGeometry& occgeo,
2106 netgen::Mesh& ngMesh,
2107 const NETGENPlugin_ngMeshInfo& initState,
2109 std::vector<const SMDS_MeshNode*>& nodeVec,
2110 SMESH_Comment& comment,
2111 SMESH_MesherHelper* quadHelper)
2113 int nbNod = ngMesh.GetNP();
2114 int nbSeg = ngMesh.GetNSeg();
2115 int nbFac = ngMesh.GetNSE();
2116 int nbVol = ngMesh.GetNE();
2118 SMESHDS_Mesh* meshDS = sMesh.GetMeshDS();
2120 // quadHelper is used for either
2121 // 1) making quadratic elements when a lower dimention mesh is loaded
2122 // to SMESH before convertion to quadratic by NETGEN
2123 // 2) sewing of quadratic elements with quadratic elements of sub-meshes
2124 if ( quadHelper && !quadHelper->GetIsQuadratic() && quadHelper->GetTLinkNodeMap().empty() )
2127 // -------------------------------------
2128 // Create and insert nodes into nodeVec
2129 // -------------------------------------
2131 nodeVec.resize( nbNod + 1 );
2132 int i, nbInitNod = initState._nbNodes;
2133 for (i = nbInitNod+1; i <= nbNod; ++i )
2135 const netgen::MeshPoint& ngPoint = ngMesh.Point(i);
2136 SMDS_MeshNode* node = NULL;
2137 TopoDS_Vertex aVert;
2138 // First, netgen creates nodes on vertices in occgeo.vmap,
2139 // so node index corresponds to vertex index
2140 // but (issue 0020776) netgen does not create nodes with equal coordinates
2141 if ( i-nbInitNod <= occgeo.vmap.Extent() )
2143 gp_Pnt p ( NGPOINT_COORDS(ngPoint) );
2144 for (int iV = i-nbInitNod; aVert.IsNull() && iV <= occgeo.vmap.Extent(); ++iV)
2146 aVert = TopoDS::Vertex( occgeo.vmap( iV ));
2147 gp_Pnt pV = BRep_Tool::Pnt( aVert );
2148 if ( p.SquareDistance( pV ) > 1e-20 )
2151 node = const_cast<SMDS_MeshNode*>( SMESH_Algo::VertexNode( aVert, meshDS ));
2154 if (!node) // node not found on vertex
2156 node = meshDS->AddNode( NGPOINT_COORDS( ngPoint ));
2157 if (!aVert.IsNull())
2158 meshDS->SetNodeOnVertex(node, aVert);
2163 // -------------------------------------------
2164 // Create mesh segments along geometric edges
2165 // -------------------------------------------
2167 int nbInitSeg = initState._nbSegments;
2168 for (i = nbInitSeg+1; i <= nbSeg; ++i )
2170 const netgen::Segment& seg = ngMesh.LineSegment(i);
2172 int pinds[3] = { seg.pnums[0], seg.pnums[1], seg.pnums[2] };
2175 for (int j=0; j < 3; ++j)
2177 int pind = pinds[j];
2178 if (pind <= 0 || !nodeVec_ACCESS(pind))
2186 int aGeomEdgeInd = seg.epgeominfo[j].edgenr;
2187 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
2188 aEdge = TopoDS::Edge(occgeo.emap(aGeomEdgeInd));
2190 param = seg.epgeominfo[j].dist;
2193 else // middle point
2195 param = param2 * 0.5;
2197 if (!aEdge.IsNull() && nodeVec_ACCESS(pind)->getshapeId() < 1)
2199 meshDS->SetNodeOnEdge(nodeVec_ACCESS(pind), aEdge, param);
2204 SMDS_MeshEdge* edge = 0;
2205 if (nbp == 2) // second order ?
2207 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1])))
2209 if ( quadHelper ) // final mesh must be quadratic
2210 edge = quadHelper->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
2212 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
2216 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
2217 nodeVec_ACCESS(pinds[2])))
2219 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
2220 nodeVec_ACCESS(pinds[2]));
2224 if ( comment.empty() ) comment << "Cannot create a mesh edge";
2225 MESSAGE("Cannot create a mesh edge");
2226 nbSeg = nbFac = nbVol = 0;
2229 if ( !aEdge.IsNull() && edge->getshapeId() < 1 )
2230 meshDS->SetMeshElementOnShape(edge, aEdge);
2232 else if ( comment.empty() )
2234 comment << "Invalid netgen segment #" << i;
2238 // ----------------------------------------
2239 // Create mesh faces along geometric faces
2240 // ----------------------------------------
2242 int nbInitFac = initState._nbFaces;
2243 int quadFaceID = ngMesh.GetNFD() + 1;
2244 if ( nbInitFac < nbFac )
2245 // add a faces descriptor to exclude qudrangle elements generated by NETGEN
2246 // from computation of 3D mesh
2247 ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(quadFaceID, /*solid1=*/0, /*solid2=*/0, 0));
2249 vector<const SMDS_MeshNode*> nodes;
2250 for (i = nbInitFac+1; i <= nbFac; ++i )
2252 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
2253 const int aGeomFaceInd = elem.GetIndex();
2255 if (aGeomFaceInd > 0 && aGeomFaceInd <= occgeo.fmap.Extent())
2256 aFace = TopoDS::Face(occgeo.fmap(aGeomFaceInd));
2258 for ( int j = 1; j <= elem.GetNP(); ++j )
2260 int pind = elem.PNum(j);
2261 if ( pind < 1 || pind >= (int) nodeVec.size() )
2263 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind))
2265 nodes.push_back( node );
2266 if (!aFace.IsNull() && node->getshapeId() < 1)
2268 const netgen::PointGeomInfo& pgi = elem.GeomInfoPi(j);
2269 meshDS->SetNodeOnFace(node, aFace, pgi.u, pgi.v);
2273 if ((int) nodes.size() != elem.GetNP() )
2275 if ( comment.empty() )
2276 comment << "Invalid netgen 2d element #" << i;
2277 continue; // bad node ids
2279 SMDS_MeshFace* face = NULL;
2280 switch (elem.GetType())
2283 if ( quadHelper ) // final mesh must be quadratic
2284 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2]);
2286 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2]);
2289 if ( quadHelper ) // final mesh must be quadratic
2290 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2292 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2293 // exclude qudrangle elements from computation of 3D mesh
2294 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2297 nodes[5] = mediumNode( nodes[0],nodes[1],nodes[5], quadHelper );
2298 nodes[3] = mediumNode( nodes[1],nodes[2],nodes[3], quadHelper );
2299 nodes[4] = mediumNode( nodes[2],nodes[0],nodes[4], quadHelper );
2300 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[5],nodes[3],nodes[4]);
2303 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2304 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2305 nodes[5] = mediumNode( nodes[2],nodes[3],nodes[5], quadHelper );
2306 nodes[6] = mediumNode( nodes[3],nodes[0],nodes[6], quadHelper );
2307 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3],
2308 nodes[4],nodes[7],nodes[5],nodes[6]);
2309 // exclude qudrangle elements from computation of 3D mesh
2310 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2313 MESSAGE("NETGEN created a face of unexpected type, ignoring");
2318 if ( comment.empty() ) comment << "Cannot create a mesh face";
2319 MESSAGE("Cannot create a mesh face");
2320 nbSeg = nbFac = nbVol = 0;
2323 if ( !aFace.IsNull() )
2324 meshDS->SetMeshElementOnShape( face, aFace );
2327 // ------------------
2328 // Create tetrahedra
2329 // ------------------
2331 for ( i = 1; i <= nbVol; ++i )
2333 const netgen::Element& elem = ngMesh.VolumeElement(i);
2334 int aSolidInd = elem.GetIndex();
2335 TopoDS_Solid aSolid;
2336 if ( aSolidInd > 0 && aSolidInd <= occgeo.somap.Extent() )
2337 aSolid = TopoDS::Solid(occgeo.somap(aSolidInd));
2339 for ( int j = 1; j <= elem.GetNP(); ++j )
2341 int pind = elem.PNum(j);
2342 if ( pind < 1 || pind >= (int)nodeVec.size() )
2344 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind) )
2346 nodes.push_back(node);
2347 if ( !aSolid.IsNull() && node->getshapeId() < 1 )
2348 meshDS->SetNodeInVolume(node, aSolid);
2351 if ((int) nodes.size() != elem.GetNP() )
2353 if ( comment.empty() )
2354 comment << "Invalid netgen 3d element #" << i;
2357 SMDS_MeshVolume* vol = NULL;
2358 switch ( elem.GetType() )
2361 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3]);
2364 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2365 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2366 nodes[5] = mediumNode( nodes[2],nodes[0],nodes[5], quadHelper );
2367 nodes[6] = mediumNode( nodes[0],nodes[3],nodes[6], quadHelper );
2368 nodes[8] = mediumNode( nodes[1],nodes[3],nodes[8], quadHelper );
2369 nodes[9] = mediumNode( nodes[2],nodes[3],nodes[9], quadHelper );
2370 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3],
2371 nodes[4],nodes[7],nodes[5],nodes[6],nodes[8],nodes[9]);
2374 MESSAGE("NETGEN created a volume of unexpected type, ignoring");
2379 if ( comment.empty() ) comment << "Cannot create a mesh volume";
2380 MESSAGE("Cannot create a mesh volume");
2381 nbSeg = nbFac = nbVol = 0;
2384 if (!aSolid.IsNull())
2385 meshDS->SetMeshElementOnShape(vol, aSolid);
2387 return comment.empty() ? 0 : 1;
2392 //================================================================================
2394 * \brief Convert error into text
2396 //================================================================================
2398 std::string text(int err)
2403 SMESH_Comment("Error in netgen::OCCGenerateMesh() at ") << netgen::multithread.task;
2406 //================================================================================
2408 * \brief Convert exception into text
2410 //================================================================================
2412 std::string text(Standard_Failure& ex)
2414 SMESH_Comment str("Exception in netgen::OCCGenerateMesh()");
2415 str << " at " << netgen::multithread.task
2416 << ": " << ex.DynamicType()->Name();
2417 if ( ex.GetMessageString() && strlen( ex.GetMessageString() ))
2418 str << ": " << ex.GetMessageString();
2421 //================================================================================
2423 * \brief Convert exception into text
2425 //================================================================================
2427 std::string text(netgen::NgException& ex)
2429 SMESH_Comment str("NgException");
2430 if ( strlen( netgen::multithread.task ) > 0 )
2431 str << " at " << netgen::multithread.task;
2432 str << ": " << ex.What();
2436 //================================================================================
2438 * \brief Looks for triangles lying on a SOLID
2440 //================================================================================
2442 bool hasBadElemOnSolid( const list<const SMDS_MeshElement*>& elems,
2443 SMESH_subMesh* solidSM )
2445 TopTools_IndexedMapOfShape solidSubs;
2446 TopExp::MapShapes( solidSM->GetSubShape(), solidSubs );
2447 SMESHDS_Mesh* mesh = solidSM->GetFather()->GetMeshDS();
2449 list<const SMDS_MeshElement*>::const_iterator e = elems.begin();
2450 for ( ; e != elems.end(); ++e )
2452 const SMDS_MeshElement* elem = *e;
2453 // if ( elem->GetType() != SMDSAbs_Face ) -- 23047
2455 int nbNodesOnSolid = 0, nbNodes = elem->NbNodes();
2456 SMDS_NodeIteratorPtr nIt = elem->nodeIterator();
2457 while ( nIt->more() )
2459 const SMDS_MeshNode* n = nIt->next();
2460 const TopoDS_Shape& s = mesh->IndexToShape( n->getshapeId() );
2461 nbNodesOnSolid += ( !s.IsNull() && solidSubs.Contains( s ));
2462 if ( nbNodesOnSolid > 2 ||
2463 nbNodesOnSolid == nbNodes)
2470 const double edgeMeshingTime = 0.001;
2471 const double faceMeshingTime = 0.019;
2472 const double edgeFaceMeshingTime = edgeMeshingTime + faceMeshingTime;
2473 const double faceOptimizTime = 0.06;
2474 const double voluMeshingTime = 0.15;
2475 const double volOptimizeTime = 0.77;
2478 //=============================================================================
2480 * Here we are going to use the NETGEN mesher
2482 //=============================================================================
2484 bool NETGENPlugin_Mesher::Compute()
2486 NETGENPlugin_NetgenLibWrapper ngLib;
2488 netgen::MeshingParameters& mparams = netgen::mparam;
2489 MESSAGE("Compute with:\n"
2490 " max size = " << mparams.maxh << "\n"
2491 " segments per edge = " << mparams.segmentsperedge);
2493 " growth rate = " << mparams.grading << "\n"
2494 " elements per radius = " << mparams.curvaturesafety << "\n"
2495 " second order = " << mparams.secondorder << "\n"
2496 " quad allowed = " << mparams.quad << "\n"
2497 " surface curvature = " << mparams.uselocalh << "\n"
2498 " fuse edges = " << netgen::merge_solids);
2500 SMESH_ComputeErrorPtr error = SMESH_ComputeError::New();
2501 SMESH_MesherHelper quadHelper( *_mesh );
2502 quadHelper.SetIsQuadratic( mparams.secondorder );
2504 static string debugFile = "/tmp/ngMesh.py"; /* to call toPython( _ngMesh, debugFile )
2505 while debugging netgen */
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 err = 0; //- MESHCONST_ANALYSE isn't so important step
2591 ngLib.setMesh(( Ng_Mesh*) _ngMesh );
2593 _ngMesh->ClearFaceDescriptors(); // we make descriptors our-self
2597 // Pass 1D simple parameters to NETGEN
2598 // --------------------------------
2599 int nbSeg = _simpleHyp->GetNumberOfSegments();
2600 double segSize = _simpleHyp->GetLocalLength();
2601 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
2603 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
2605 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
2606 setLocalSize( e, segSize, *_ngMesh );
2609 else // if ( ! _simpleHyp )
2611 // Local size on shapes
2612 SetLocalSize( occgeo, *_ngMesh );
2615 // Precompute internal edges (issue 0020676) in order to
2616 // add mesh on them correctly (twice) to netgen mesh
2617 if ( !err && internals.hasInternalEdges() )
2619 // load internal shapes into OCCGeometry
2620 netgen::OCCGeometry intOccgeo;
2621 internals.getInternalEdges( intOccgeo.fmap, intOccgeo.emap, intOccgeo.vmap, meshedSM );
2622 intOccgeo.boundingbox = occgeo.boundingbox;
2623 intOccgeo.shape = occgeo.shape;
2624 intOccgeo.face_maxh.SetSize(intOccgeo.fmap.Extent());
2625 intOccgeo.face_maxh = netgen::mparam.maxh;
2626 netgen::Mesh *tmpNgMesh = NULL;
2630 // compute local H on internal shapes in the main mesh
2631 //OCCSetLocalMeshSize(intOccgeo, *_ngMesh); it deletes _ngMesh->localH
2633 // let netgen create a temporary mesh
2635 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2637 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2639 if(netgen::multithread.terminate)
2642 // copy LocalH from the main to temporary mesh
2643 initState.transferLocalH( _ngMesh, tmpNgMesh );
2645 // compute mesh on internal edges
2646 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2648 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2650 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2652 comment << text(err);
2654 catch (Standard_Failure& ex)
2656 comment << text(ex);
2659 initState.restoreLocalH( tmpNgMesh );
2661 // fill SMESH by netgen mesh
2662 vector< const SMDS_MeshNode* > tmpNodeVec;
2663 FillSMesh( intOccgeo, *tmpNgMesh, initState, *_mesh, tmpNodeVec, comment );
2664 err = ( err || !comment.empty() );
2666 nglib::Ng_DeleteMesh((nglib::Ng_Mesh*)tmpNgMesh);
2669 // Fill _ngMesh with nodes and segments of computed submeshes
2672 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_0D ]) &&
2673 FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_1D ], &quadHelper));
2675 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2680 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2685 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2687 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2689 if(netgen::multithread.terminate)
2692 comment << text(err);
2694 catch (Standard_Failure& ex)
2696 comment << text(ex);
2701 _ticTime = ( doneTime += edgeMeshingTime ) / _totalTime / _progressTic;
2703 mparams.uselocalh = true; // restore as it is used at surface optimization
2705 // ---------------------
2706 // compute surface mesh
2707 // ---------------------
2710 // Pass 2D simple parameters to NETGEN
2712 if ( double area = _simpleHyp->GetMaxElementArea() ) {
2714 mparams.maxh = sqrt(2. * area/sqrt(3.0));
2715 mparams.grading = 0.4; // moderate size growth
2718 // length from edges
2719 if ( _ngMesh->GetNSeg() ) {
2720 double edgeLength = 0;
2721 TopTools_MapOfShape visitedEdges;
2722 for ( TopExp_Explorer exp( _shape, TopAbs_EDGE ); exp.More(); exp.Next() )
2723 if( visitedEdges.Add(exp.Current()) )
2724 edgeLength += SMESH_Algo::EdgeLength( TopoDS::Edge( exp.Current() ));
2725 // we have to multiply length by 2 since for each TopoDS_Edge there
2726 // are double set of NETGEN edges, in other words, we have to
2727 // divide _ngMesh->GetNSeg() by 2.
2728 mparams.maxh = 2*edgeLength / _ngMesh->GetNSeg();
2731 mparams.maxh = 1000;
2733 mparams.grading = 0.2; // slow size growth
2735 mparams.quad = _simpleHyp->GetAllowQuadrangles();
2736 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2737 _ngMesh->SetGlobalH (mparams.maxh);
2738 netgen::Box<3> bb = occgeo.GetBoundingBox();
2739 bb.Increase (bb.Diam()/20);
2740 _ngMesh->SetLocalH (bb.PMin(), bb.PMax(), mparams.grading);
2743 // Care of vertices internal in faces (issue 0020676)
2744 if ( internals.hasInternalVertexInFace() )
2746 // store computed segments in SMESH in order not to create SMESH
2747 // edges for ng segments added by AddIntVerticesInFaces()
2748 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2749 // add segments to faces with internal vertices
2750 AddIntVerticesInFaces( occgeo, *_ngMesh, nodeVec, internals );
2751 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2754 // Build viscous layers
2755 if ( _isViscousLayers2D ||
2756 StdMeshers_ViscousLayers2D::HasProxyMesh( TopoDS::Face( occgeo.fmap(1) ), *_mesh ))
2758 if ( !internals.hasInternalVertexInFace() ) {
2759 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2760 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2762 SMESH_ProxyMesh::Ptr viscousMesh;
2763 SMESH_MesherHelper helper( *_mesh );
2764 for ( int faceID = 1; faceID <= occgeo.fmap.Extent(); ++faceID )
2766 const TopoDS_Face& F = TopoDS::Face( occgeo.fmap( faceID ));
2767 viscousMesh = StdMeshers_ViscousLayers2D::Compute( *_mesh, F );
2770 if ( viscousMesh->NbProxySubMeshes() == 0 )
2772 // exclude from computation ng segments built on EDGEs of F
2773 for (int i = 1; i <= _ngMesh->GetNSeg(); i++)
2775 netgen::Segment & seg = _ngMesh->LineSegment(i);
2776 if (seg.si == faceID)
2779 // add new segments to _ngMesh instead of excluded ones
2780 helper.SetSubShape( F );
2782 StdMeshers_FaceSide::GetFaceWires( F, *_mesh, /*skipMediumNodes=*/true,
2783 error, viscousMesh );
2784 error = AddSegmentsToMesh( *_ngMesh, occgeo, wires, helper, nodeVec );
2786 if ( !error ) error = SMESH_ComputeError::New();
2788 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2791 // Let netgen compute 2D mesh
2792 startWith = netgen::MESHCONST_MESHSURFACE;
2793 endWith = _optimize ? netgen::MESHCONST_OPTSURFACE : netgen::MESHCONST_MESHSURFACE;
2798 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2800 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2802 if(netgen::multithread.terminate)
2805 comment << text (err);
2807 catch (Standard_Failure& ex)
2809 comment << text(ex);
2810 //err = 1; -- try to make volumes anyway
2812 catch (netgen::NgException exc)
2814 comment << text(exc);
2815 //err = 1; -- try to make volumes anyway
2820 doneTime += faceMeshingTime + ( _optimize ? faceOptimizTime : 0 );
2821 _ticTime = doneTime / _totalTime / _progressTic;
2823 // ---------------------
2824 // generate volume mesh
2825 // ---------------------
2826 // Fill _ngMesh with nodes and faces of computed 2D submeshes
2827 if ( !err && _isVolume && ( !meshedSM[ MeshDim_2D ].empty() || mparams.quad ))
2829 // load SMESH with computed segments and faces
2830 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2832 // compute pyramids on quadrangles
2833 SMESH_ProxyMesh::Ptr proxyMesh;
2834 if ( _mesh->NbQuadrangles() > 0 )
2835 for ( int iS = 1; iS <= occgeo.somap.Extent(); ++iS )
2837 StdMeshers_QuadToTriaAdaptor* Adaptor = new StdMeshers_QuadToTriaAdaptor;
2838 proxyMesh.reset( Adaptor );
2840 int nbPyrams = _mesh->NbPyramids();
2841 Adaptor->Compute( *_mesh, occgeo.somap(iS) );
2842 if ( nbPyrams != _mesh->NbPyramids() )
2844 list< SMESH_subMesh* > quadFaceSM;
2845 for (TopExp_Explorer face(occgeo.somap(iS), TopAbs_FACE); face.More(); face.Next())
2846 if ( Adaptor->GetProxySubMesh( face.Current() ))
2848 quadFaceSM.push_back( _mesh->GetSubMesh( face.Current() ));
2849 meshedSM[ MeshDim_2D ].remove( quadFaceSM.back() );
2851 FillNgMesh(occgeo, *_ngMesh, nodeVec, quadFaceSM, &quadHelper, proxyMesh);
2854 // fill _ngMesh with faces of sub-meshes
2855 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_2D ], &quadHelper));
2856 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2857 //toPython( _ngMesh, "/tmp/ngPython.py");
2859 if (!err && _isVolume)
2861 // Pass 3D simple parameters to NETGEN
2862 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
2863 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
2865 if ( double vol = simple3d->GetMaxElementVolume() ) {
2867 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
2868 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2871 // length from faces
2872 mparams.maxh = _ngMesh->AverageH();
2874 _ngMesh->SetGlobalH (mparams.maxh);
2875 mparams.grading = 0.4;
2877 _ngMesh->CalcLocalH(mparams.grading);
2879 _ngMesh->CalcLocalH();
2882 // Care of vertices internal in solids and internal faces (issue 0020676)
2883 if ( internals.hasInternalVertexInSolid() || internals.hasInternalFaces() )
2885 // store computed faces in SMESH in order not to create SMESH
2886 // faces for ng faces added here
2887 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2888 // add ng faces to solids with internal vertices
2889 AddIntVerticesInSolids( occgeo, *_ngMesh, nodeVec, internals );
2890 // duplicate mesh faces on internal faces
2891 FixIntFaces( occgeo, *_ngMesh, internals );
2892 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2894 // Let netgen compute 3D mesh
2895 startWith = endWith = netgen::MESHCONST_MESHVOLUME;
2900 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2902 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2904 if(netgen::multithread.terminate)
2907 if ( comment.empty() ) // do not overwrite a previos error
2908 comment << text(err);
2910 catch (Standard_Failure& ex)
2912 if ( comment.empty() ) // do not overwrite a previos error
2913 comment << text(ex);
2916 catch (netgen::NgException exc)
2918 if ( comment.empty() ) // do not overwrite a previos error
2919 comment << text(exc);
2922 _ticTime = ( doneTime += voluMeshingTime ) / _totalTime / _progressTic;
2924 // Let netgen optimize 3D mesh
2925 if ( !err && _optimize )
2927 startWith = endWith = netgen::MESHCONST_OPTVOLUME;
2932 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2934 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2936 if(netgen::multithread.terminate)
2939 if ( comment.empty() ) // do not overwrite a previos error
2940 comment << text(err);
2942 catch (Standard_Failure& ex)
2944 if ( comment.empty() ) // do not overwrite a previos error
2945 comment << text(ex);
2947 catch (netgen::NgException exc)
2949 if ( comment.empty() ) // do not overwrite a previos error
2950 comment << text(exc);
2954 if (!err && mparams.secondorder > 0)
2959 if ( !meshedSM[ MeshDim_1D ].empty() )
2961 // remove segments not attached to geometry (IPAL0052479)
2962 for (int i = 1; i <= _ngMesh->GetNSeg(); ++i)
2964 const netgen::Segment & seg = _ngMesh->LineSegment (i);
2965 if ( seg.epgeominfo[ 0 ].edgenr == 0 )
2966 _ngMesh->DeleteSegment( i );
2968 _ngMesh->Compress();
2970 // convert to quadratic
2971 netgen::OCCRefinementSurfaces ref (occgeo);
2972 ref.MakeSecondOrder (*_ngMesh);
2974 // care of elements already loaded to SMESH
2975 // if ( initState._nbSegments > 0 )
2976 // makeQuadratic( occgeo.emap, _mesh );
2977 // if ( initState._nbFaces > 0 )
2978 // makeQuadratic( occgeo.fmap, _mesh );
2980 catch (Standard_Failure& ex)
2982 if ( comment.empty() ) // do not overwrite a previos error
2983 comment << "Exception in netgen at passing to 2nd order ";
2985 catch (netgen::NgException exc)
2987 if ( comment.empty() ) // do not overwrite a previos error
2988 comment << exc.What();
2993 _ticTime = 0.98 / _progressTic;
2995 //int nbNod = _ngMesh->GetNP();
2996 //int nbSeg = _ngMesh->GetNSeg();
2997 int nbFac = _ngMesh->GetNSE();
2998 int nbVol = _ngMesh->GetNE();
2999 bool isOK = ( !err && (_isVolume ? (nbVol > 0) : (nbFac > 0)) );
3001 // Feed back the SMESHDS with the generated Nodes and Elements
3002 if ( true /*isOK*/ ) // get whatever built
3004 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
3006 if ( quadHelper.GetIsQuadratic() ) // remove free nodes
3007 for ( size_t i = 0; i < nodeVec.size(); ++i )
3008 if ( nodeVec[i] && nodeVec[i]->NbInverseElements() == 0 )
3009 _mesh->GetMeshDS()->RemoveFreeNode( nodeVec[i], 0, /*fromGroups=*/false );
3011 SMESH_ComputeErrorPtr readErr = ReadErrors(nodeVec);
3012 if ( readErr && !readErr->myBadElements.empty() )
3015 if ( !comment.empty() && !readErr->myComment.empty() ) comment += "\n";
3016 comment += readErr->myComment;
3018 if ( error->IsOK() && ( !isOK || comment.size() > 0 ))
3019 error->myName = COMPERR_ALGO_FAILED;
3020 if ( !comment.empty() )
3021 error->myComment = comment;
3023 // SetIsAlwaysComputed( true ) to empty sub-meshes, which
3024 // appear if the geometry contains coincident sub-shape due
3025 // to bool merge_solids = 1; in netgen/libsrc/occ/occgenmesh.cpp
3026 const int nbMaps = 2;
3027 const TopTools_IndexedMapOfShape* geoMaps[nbMaps] =
3028 { & occgeo.vmap, & occgeo.emap/*, & occgeo.fmap*/ };
3029 for ( int iMap = 0; iMap < nbMaps; ++iMap )
3030 for (int i = 1; i <= geoMaps[iMap]->Extent(); i++)
3031 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( geoMaps[iMap]->FindKey(i)))
3032 if ( !sm->IsMeshComputed() )
3033 sm->SetIsAlwaysComputed( true );
3035 // set bad compute error to subshapes of all failed sub-shapes
3036 if ( !error->IsOK() )
3038 bool pb2D = false, pb3D = false;
3039 for (int i = 1; i <= occgeo.fmap.Extent(); i++) {
3040 int status = occgeo.facemeshstatus[i-1];
3041 if (status == netgen::FACE_MESHED_OK ) continue;
3042 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.fmap( i ))) {
3043 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
3044 if ( !smError || smError->IsOK() ) {
3045 if ( status == netgen::FACE_FAILED )
3046 smError.reset( new SMESH_ComputeError( *error ));
3048 smError.reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED, "Ignored" ));
3049 if ( SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
3050 smError->myName = COMPERR_WARNING;
3052 pb2D = pb2D || smError->IsKO();
3055 if ( !pb2D ) // all faces are OK
3056 for (int i = 1; i <= occgeo.somap.Extent(); i++)
3057 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.somap( i )))
3059 bool smComputed = nbVol && !sm->IsEmpty();
3060 if ( smComputed && internals.hasInternalVertexInSolid( sm->GetId() ))
3062 int nbIntV = internals.getSolidsWithVertices().find( sm->GetId() )->second.size();
3063 SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
3064 smComputed = ( smDS->NbElements() > 0 || smDS->NbNodes() > nbIntV );
3066 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
3067 if ( !smComputed && ( !smError || smError->IsOK() ))
3069 smError.reset( new SMESH_ComputeError( *error ));
3070 if ( nbVol && SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
3072 smError->myName = COMPERR_WARNING;
3074 else if ( !smError->myBadElements.empty() ) // bad surface mesh
3076 if ( !hasBadElemOnSolid( smError->myBadElements, sm ))
3080 pb3D = pb3D || ( smError && smError->IsKO() );
3082 if ( !pb2D && !pb3D )
3083 err = 0; // no fatal errors, only warnings
3086 ngLib._isComputeOk = !err;
3091 //=============================================================================
3095 //=============================================================================
3096 bool NETGENPlugin_Mesher::Evaluate(MapShapeNbElems& aResMap)
3098 netgen::MeshingParameters& mparams = netgen::mparam;
3101 // -------------------------
3102 // Prepare OCC geometry
3103 // -------------------------
3104 netgen::OCCGeometry occgeo;
3105 list< SMESH_subMesh* > meshedSM[4]; // for 0-3 dimensions
3106 NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
3107 PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
3109 bool tooManyElems = false;
3110 const int hugeNb = std::numeric_limits<int>::max() / 100;
3115 // pass 1D simple parameters to NETGEN
3118 // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
3119 mparams.uselocalh = false;
3120 mparams.grading = 0.8; // not limitited size growth
3122 if ( _simpleHyp->GetNumberOfSegments() )
3124 mparams.maxh = occgeo.boundingbox.Diam();
3127 mparams.maxh = _simpleHyp->GetLocalLength();
3130 if ( mparams.maxh == 0.0 )
3131 mparams.maxh = occgeo.boundingbox.Diam();
3132 if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
3133 mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
3135 // let netgen create _ngMesh and calculate element size on not meshed shapes
3136 NETGENPlugin_NetgenLibWrapper ngLib;
3137 netgen::Mesh *ngMesh = NULL;
3141 int startWith = netgen::MESHCONST_ANALYSE;
3142 int endWith = netgen::MESHCONST_MESHEDGES;
3144 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, mparams, startWith, endWith);
3146 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, startWith, endWith, optstr);
3149 if(netgen::multithread.terminate)
3152 ngLib.setMesh(( Ng_Mesh*) ngMesh );
3154 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( _shape ))
3155 sm->GetComputeError().reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED ));
3160 // Pass 1D simple parameters to NETGEN
3161 // --------------------------------
3162 int nbSeg = _simpleHyp->GetNumberOfSegments();
3163 double segSize = _simpleHyp->GetLocalLength();
3164 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
3166 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
3168 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
3169 setLocalSize( e, segSize, *ngMesh );
3172 else // if ( ! _simpleHyp )
3174 // Local size on shapes
3175 SetLocalSize( occgeo, *ngMesh );
3177 // calculate total nb of segments and length of edges
3178 double fullLen = 0.0;
3180 int entity = mparams.secondorder > 0 ? SMDSEntity_Quad_Edge : SMDSEntity_Edge;
3181 TopTools_DataMapOfShapeInteger Edge2NbSeg;
3182 for (TopExp_Explorer exp(_shape, TopAbs_EDGE); exp.More(); exp.Next())
3184 TopoDS_Edge E = TopoDS::Edge( exp.Current() );
3185 if( !Edge2NbSeg.Bind(E,0) )
3188 double aLen = SMESH_Algo::EdgeLength(E);
3191 vector<int>& aVec = aResMap[_mesh->GetSubMesh(E)];
3193 aVec.resize( SMDSEntity_Last, 0);
3195 fullNbSeg += aVec[ entity ];
3198 // store nb of segments computed by Netgen
3199 NCollection_Map<Link> linkMap;
3200 for (int i = 1; i <= ngMesh->GetNSeg(); ++i )
3202 const netgen::Segment& seg = ngMesh->LineSegment(i);
3203 Link link(seg[0], seg[1]);
3204 if ( !linkMap.Add( link )) continue;
3205 int aGeomEdgeInd = seg.epgeominfo[0].edgenr;
3206 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
3208 vector<int>& aVec = aResMap[_mesh->GetSubMesh(occgeo.emap(aGeomEdgeInd))];
3212 // store nb of nodes on edges computed by Netgen
3213 TopTools_DataMapIteratorOfDataMapOfShapeInteger Edge2NbSegIt(Edge2NbSeg);
3214 for (; Edge2NbSegIt.More(); Edge2NbSegIt.Next())
3216 vector<int>& aVec = aResMap[_mesh->GetSubMesh(Edge2NbSegIt.Key())];
3217 if ( aVec[ entity ] > 1 && aVec[ SMDSEntity_Node ] == 0 )
3218 aVec[SMDSEntity_Node] = mparams.secondorder > 0 ? 2*aVec[ entity ]-1 : aVec[ entity ]-1;
3220 fullNbSeg += aVec[ entity ];
3221 Edge2NbSeg( Edge2NbSegIt.Key() ) = aVec[ entity ];
3223 if ( fullNbSeg == 0 )
3230 if ( double area = _simpleHyp->GetMaxElementArea() ) {
3232 mparams.maxh = sqrt(2. * area/sqrt(3.0));
3233 mparams.grading = 0.4; // moderate size growth
3236 // length from edges
3237 mparams.maxh = fullLen/fullNbSeg;
3238 mparams.grading = 0.2; // slow size growth
3241 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3242 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3244 for (TopExp_Explorer exp(_shape, TopAbs_FACE); exp.More(); exp.Next())
3246 TopoDS_Face F = TopoDS::Face( exp.Current() );
3247 SMESH_subMesh *sm = _mesh->GetSubMesh(F);
3249 BRepGProp::SurfaceProperties(F,G);
3250 double anArea = G.Mass();
3251 tooManyElems = tooManyElems || ( anArea/hugeNb > mparams.maxh*mparams.maxh );
3253 if ( !tooManyElems )
3255 TopTools_MapOfShape egdes;
3256 for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next())
3257 if ( egdes.Add( exp1.Current() ))
3258 nb1d += Edge2NbSeg.Find(exp1.Current());
3260 int nbFaces = tooManyElems ? hugeNb : int( 4*anArea / (mparams.maxh*mparams.maxh*sqrt(3.)));
3261 int nbNodes = tooManyElems ? hugeNb : (( nbFaces*3 - (nb1d-1)*2 ) / 6 + 1 );
3263 vector<int> aVec(SMDSEntity_Last, 0);
3264 if( mparams.secondorder > 0 ) {
3265 int nb1d_in = (nbFaces*3 - nb1d) / 2;
3266 aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
3267 aVec[SMDSEntity_Quad_Triangle] = nbFaces;
3270 aVec[SMDSEntity_Node] = Max ( nbNodes, 0 );
3271 aVec[SMDSEntity_Triangle] = nbFaces;
3273 aResMap[sm].swap(aVec);
3280 // pass 3D simple parameters to NETGEN
3281 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
3282 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
3284 if ( double vol = simple3d->GetMaxElementVolume() ) {
3286 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
3287 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3290 // using previous length from faces
3292 mparams.grading = 0.4;
3293 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3296 BRepGProp::VolumeProperties(_shape,G);
3297 double aVolume = G.Mass();
3298 double tetrVol = 0.1179*mparams.maxh*mparams.maxh*mparams.maxh;
3299 tooManyElems = tooManyElems || ( aVolume/hugeNb > tetrVol );
3300 int nbVols = tooManyElems ? hugeNb : int(aVolume/tetrVol);
3301 int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
3302 vector<int> aVec(SMDSEntity_Last, 0 );
3303 if ( tooManyElems ) // avoid FPE
3305 aVec[SMDSEntity_Node] = hugeNb;
3306 aVec[ mparams.secondorder > 0 ? SMDSEntity_Quad_Tetra : SMDSEntity_Tetra] = hugeNb;
3310 if( mparams.secondorder > 0 ) {
3311 aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
3312 aVec[SMDSEntity_Quad_Tetra] = nbVols;
3315 aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
3316 aVec[SMDSEntity_Tetra] = nbVols;
3319 SMESH_subMesh *sm = _mesh->GetSubMesh(_shape);
3320 aResMap[sm].swap(aVec);
3326 double NETGENPlugin_Mesher::GetProgress(const SMESH_Algo* holder,
3327 const int * algoProgressTic,
3328 const double * algoProgress) const
3330 ((int&) _progressTic ) = *algoProgressTic + 1;
3332 if ( !_occgeom ) return 0;
3334 double progress = -1;
3337 if ( _ticTime < 0 && netgen::multithread.task[0] == 'O'/*Optimizing surface*/ )
3339 ((double&) _ticTime ) = edgeFaceMeshingTime / _totalTime / _progressTic;
3341 else if ( !_optimize /*&& _occgeom->fmap.Extent() > 1*/ )
3343 int doneShapeIndex = -1;
3344 while ( doneShapeIndex+1 < _occgeom->facemeshstatus.Size() &&
3345 _occgeom->facemeshstatus[ doneShapeIndex+1 ])
3347 if ( doneShapeIndex+1 != _curShapeIndex )
3349 ((int&) _curShapeIndex) = doneShapeIndex+1;
3350 double doneShapeRate = _curShapeIndex / double( _occgeom->fmap.Extent() );
3351 double doneTime = edgeMeshingTime + doneShapeRate * faceMeshingTime;
3352 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3353 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3354 // << " " << doneTime / _totalTime / _progressTic << endl;
3358 else if ( !_optimize && _occgeom->somap.Extent() > 1 )
3360 int curShapeIndex = _curShapeIndex;
3361 if ( _ngMesh->GetNE() > 0 )
3363 netgen::Element el = (*_ngMesh)[netgen::ElementIndex( _ngMesh->GetNE()-1 )];
3364 curShapeIndex = el.GetIndex();
3366 if ( curShapeIndex != _curShapeIndex )
3368 ((int&) _curShapeIndex) = curShapeIndex;
3369 double doneShapeRate = _curShapeIndex / double( _occgeom->somap.Extent() );
3370 double doneTime = edgeFaceMeshingTime + doneShapeRate * voluMeshingTime;
3371 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3372 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3373 // << " " << doneTime / _totalTime / _progressTic << endl;
3378 progress = Max( *algoProgressTic * _ticTime, *algoProgress );
3383 netgen::multithread.task[0] == 'D'/*elaunay meshing*/ &&
3384 progress > voluMeshingTime )
3386 progress = voluMeshingTime;
3387 ((double&) _ticTime) = voluMeshingTime / _totalTime / _progressTic;
3389 ((int&) *algoProgressTic )++;
3390 ((double&) *algoProgress) = progress;
3392 //cout << progress << " " << *algoProgressTic << " " << netgen::multithread.task << " "<< _ticTime << endl;
3394 return Min( progress, 0.99 );
3397 //================================================================================
3399 * \brief Read mesh entities preventing successful computation from "test.out" file
3401 //================================================================================
3403 SMESH_ComputeErrorPtr
3404 NETGENPlugin_Mesher::ReadErrors(const vector<const SMDS_MeshNode* >& nodeVec)
3406 SMESH_ComputeErrorPtr err = SMESH_ComputeError::New
3407 (COMPERR_BAD_INPUT_MESH, "Some edges multiple times in surface mesh");
3408 SMESH_File file("test.out");
3410 vector<int> three1(3), three2(3);
3411 const char* badEdgeStr = " multiple times in surface mesh";
3412 const int badEdgeStrLen = strlen( badEdgeStr );
3413 const int nbNodes = nodeVec.size();
3415 while( !file.eof() )
3417 if ( strncmp( file, "Edge ", 5 ) == 0 &&
3418 file.getInts( two ) &&
3419 strncmp( file, badEdgeStr, badEdgeStrLen ) == 0 &&
3420 two[0] < nbNodes && two[1] < nbNodes )
3422 err->myBadElements.push_back( new SMDS_LinearEdge( nodeVec[ two[0]], nodeVec[ two[1]] ));
3423 file += badEdgeStrLen;
3425 else if ( strncmp( file, "Intersecting: ", 14 ) == 0 )
3428 // openelement 18 with open element 126
3432 const char* pos = file;
3433 bool ok = ( strncmp( file, "openelement ", 12 ) == 0 );
3434 ok = ok && file.getInts( two );
3435 ok = ok && file.getInts( three1 );
3436 ok = ok && file.getInts( three2 );
3437 for ( int i = 0; ok && i < 3; ++i )
3438 ok = ( three1[i] < nbNodes && nodeVec[ three1[i]]);
3439 for ( int i = 0; ok && i < 3; ++i )
3440 ok = ( three2[i] < nbNodes && nodeVec[ three2[i]]);
3443 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three1[0]],
3444 nodeVec[ three1[1]],
3445 nodeVec[ three1[2]]));
3446 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three2[0]],
3447 nodeVec[ three2[1]],
3448 nodeVec[ three2[2]]));
3449 err->myComment = "Intersecting triangles";
3463 size_t nbBadElems = err->myBadElements.size();
3464 if ( nbBadElems ) nbBadElems++; // avoid warning: variable set but not used
3470 //================================================================================
3472 * \brief Write a python script creating an equivalent SALOME mesh.
3473 * This is useful to see what mesh is passed as input for the next step of mesh
3474 * generation (of mesh of higher dimension)
3476 //================================================================================
3478 void NETGENPlugin_Mesher::toPython( const netgen::Mesh* ngMesh,
3479 const std::string& pyFile)
3481 ofstream outfile(pyFile.c_str(), ios::out);
3482 if ( !outfile ) return;
3484 outfile << "import SMESH" << endl
3485 << "from salome.smesh import smeshBuilder" << endl
3486 << "smesh = smeshBuilder.New(salome.myStudy)" << endl
3487 << "mesh = smesh.Mesh()" << endl << endl;
3489 using namespace netgen;
3491 for (pi = PointIndex::BASE;
3492 pi < ngMesh->GetNP()+PointIndex::BASE; pi++)
3494 outfile << "mesh.AddNode( ";
3495 outfile << (*ngMesh)[pi](0) << ", ";
3496 outfile << (*ngMesh)[pi](1) << ", ";
3497 outfile << (*ngMesh)[pi](2) << ") ## "<< pi << endl;
3500 int nbDom = ngMesh->GetNDomains();
3501 for ( int i = 0; i < nbDom; ++i )
3502 outfile<< "grp" << i+1 << " = mesh.CreateEmptyGroup( SMESH.FACE, 'domain"<< i+1 << "')"<< endl;
3504 SurfaceElementIndex sei;
3505 for (sei = 0; sei < ngMesh->GetNSE(); sei++)
3507 outfile << "mesh.AddFace([ ";
3508 Element2d sel = (*ngMesh)[sei];
3509 for (int j = 0; j < sel.GetNP(); j++)
3510 outfile << sel[j] << ( j+1 < sel.GetNP() ? ", " : " ])");
3511 if ( sel.IsDeleted() ) outfile << " ## IsDeleted ";
3514 if ((*ngMesh)[sei].GetIndex())
3516 if ( int dom1 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainIn())
3517 outfile << "grp"<< dom1 <<".Add([ " << (int)sei+1 << " ])" << endl;
3518 if ( int dom2 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainOut())
3519 outfile << "grp"<< dom2 <<".Add([ " << (int)sei+1 << " ])" << endl;
3523 for (ElementIndex ei = 0; ei < ngMesh->GetNE(); ei++)
3525 Element el = (*ngMesh)[ei];
3526 outfile << "mesh.AddVolume([ ";
3527 for (int j = 0; j < el.GetNP(); j++)
3528 outfile << el[j] << ( j+1 < el.GetNP() ? ", " : " ])");
3532 for (int i = 1; i <= ngMesh->GetNSeg(); i++)
3534 const Segment & seg = ngMesh->LineSegment (i);
3535 outfile << "mesh.AddEdge([ "
3537 << seg[1] << " ])" << endl;
3539 cout << "Write " << pyFile << endl;
3542 //================================================================================
3544 * \brief Constructor of NETGENPlugin_ngMeshInfo
3546 //================================================================================
3548 NETGENPlugin_ngMeshInfo::NETGENPlugin_ngMeshInfo( netgen::Mesh* ngMesh):
3553 _nbNodes = ngMesh->GetNP();
3554 _nbSegments = ngMesh->GetNSeg();
3555 _nbFaces = ngMesh->GetNSE();
3556 _nbVolumes = ngMesh->GetNE();
3560 _nbNodes = _nbSegments = _nbFaces = _nbVolumes = 0;
3564 //================================================================================
3566 * \brief Copy LocalH member from one netgen mesh to another
3568 //================================================================================
3570 void NETGENPlugin_ngMeshInfo::transferLocalH( netgen::Mesh* fromMesh,
3571 netgen::Mesh* toMesh )
3573 if ( !fromMesh->LocalHFunctionGenerated() ) return;
3574 if ( !toMesh->LocalHFunctionGenerated() )
3576 toMesh->CalcLocalH(netgen::mparam.grading);
3578 toMesh->CalcLocalH();
3581 const size_t size = sizeof( netgen::LocalH );
3582 _copyOfLocalH = new char[ size ];
3583 memcpy( (void*)_copyOfLocalH, (void*)&toMesh->LocalHFunction(), size );
3584 memcpy( (void*)&toMesh->LocalHFunction(), (void*)&fromMesh->LocalHFunction(), size );
3587 //================================================================================
3589 * \brief Restore LocalH member of a netgen mesh
3591 //================================================================================
3593 void NETGENPlugin_ngMeshInfo::restoreLocalH( netgen::Mesh* toMesh )
3595 if ( _copyOfLocalH )
3597 const size_t size = sizeof( netgen::LocalH );
3598 memcpy( (void*)&toMesh->LocalHFunction(), (void*)_copyOfLocalH, size );
3599 delete [] _copyOfLocalH;
3604 //================================================================================
3606 * \brief Find "internal" sub-shapes
3608 //================================================================================
3610 NETGENPlugin_Internals::NETGENPlugin_Internals( SMESH_Mesh& mesh,
3611 const TopoDS_Shape& shape,
3613 : _mesh( mesh ), _is3D( is3D )
3615 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3617 TopExp_Explorer f,e;
3618 for ( f.Init( shape, TopAbs_FACE ); f.More(); f.Next() )
3620 int faceID = meshDS->ShapeToIndex( f.Current() );
3622 // find not computed internal edges
3624 for ( e.Init( f.Current().Oriented(TopAbs_FORWARD), TopAbs_EDGE ); e.More(); e.Next() )
3625 if ( e.Current().Orientation() == TopAbs_INTERNAL )
3627 SMESH_subMesh* eSM = mesh.GetSubMesh( e.Current() );
3628 if ( eSM->IsEmpty() )
3630 _e2face.insert( make_pair( eSM->GetId(), faceID ));
3631 for ( TopoDS_Iterator v(e.Current()); v.More(); v.Next() )
3632 _e2face.insert( make_pair( meshDS->ShapeToIndex( v.Value() ), faceID ));
3636 // find internal vertices in a face
3637 set<int> intVV; // issue 0020850 where same vertex is twice in a face
3638 for ( TopoDS_Iterator fSub( f.Current() ); fSub.More(); fSub.Next())
3639 if ( fSub.Value().ShapeType() == TopAbs_VERTEX )
3641 int vID = meshDS->ShapeToIndex( fSub.Value() );
3642 if ( intVV.insert( vID ).second )
3643 _f2v[ faceID ].push_back( vID );
3648 // find internal faces and their subshapes where nodes are to be doubled
3649 // to make a crack with non-sewed borders
3651 if ( f.Current().Orientation() == TopAbs_INTERNAL )
3653 _intShapes.insert( meshDS->ShapeToIndex( f.Current() ));
3656 list< TopoDS_Shape > edges;
3657 for ( e.Init( f.Current(), TopAbs_EDGE ); e.More(); e.Next())
3658 if ( SMESH_MesherHelper::NbAncestors( e.Current(), mesh, TopAbs_FACE ) > 1 )
3660 _intShapes.insert( meshDS->ShapeToIndex( e.Current() ));
3661 edges.push_back( e.Current() );
3662 // find border faces
3663 PShapeIteratorPtr fIt =
3664 SMESH_MesherHelper::GetAncestors( edges.back(),mesh,TopAbs_FACE );
3665 while ( const TopoDS_Shape* pFace = fIt->next() )
3666 if ( !pFace->IsSame( f.Current() ))
3667 _borderFaces.insert( meshDS->ShapeToIndex( *pFace ));
3670 // we consider vertex internal if it is shared by more than one internal edge
3671 list< TopoDS_Shape >::iterator edge = edges.begin();
3672 for ( ; edge != edges.end(); ++edge )
3673 for ( TopoDS_Iterator v( *edge ); v.More(); v.Next() )
3675 set<int> internalEdges;
3676 PShapeIteratorPtr eIt =
3677 SMESH_MesherHelper::GetAncestors( v.Value(),mesh,TopAbs_EDGE );
3678 while ( const TopoDS_Shape* pEdge = eIt->next() )
3680 int edgeID = meshDS->ShapeToIndex( *pEdge );
3681 if ( isInternalShape( edgeID ))
3682 internalEdges.insert( edgeID );
3684 if ( internalEdges.size() > 1 )
3685 _intShapes.insert( meshDS->ShapeToIndex( v.Value() ));
3689 } // loop on geom faces
3691 // find vertices internal in solids
3694 for ( TopExp_Explorer so(shape, TopAbs_SOLID); so.More(); so.Next())
3696 int soID = meshDS->ShapeToIndex( so.Current() );
3697 for ( TopoDS_Iterator soSub( so.Current() ); soSub.More(); soSub.Next())
3698 if ( soSub.Value().ShapeType() == TopAbs_VERTEX )
3699 _s2v[ soID ].push_back( meshDS->ShapeToIndex( soSub.Value() ));
3704 //================================================================================
3706 * \brief Find mesh faces on non-internal geom faces sharing internal edge
3707 * some nodes of which are to be doubled to make the second border of the "crack"
3709 //================================================================================
3711 void NETGENPlugin_Internals::findBorderElements( TIDSortedElemSet & borderElems )
3713 if ( _intShapes.empty() ) return;
3715 SMESH_Mesh& mesh = const_cast<SMESH_Mesh&>(_mesh);
3716 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3718 // loop on internal geom edges
3719 set<int>::const_iterator intShapeId = _intShapes.begin();
3720 for ( ; intShapeId != _intShapes.end(); ++intShapeId )
3722 const TopoDS_Shape& s = meshDS->IndexToShape( *intShapeId );
3723 if ( s.ShapeType() != TopAbs_EDGE ) continue;
3725 // get internal and non-internal geom faces sharing the internal edge <s>
3727 set<int>::iterator bordFace = _borderFaces.end();
3728 PShapeIteratorPtr faces = SMESH_MesherHelper::GetAncestors( s, _mesh, TopAbs_FACE );
3729 while ( const TopoDS_Shape* pFace = faces->next() )
3731 int faceID = meshDS->ShapeToIndex( *pFace );
3732 if ( isInternalShape( faceID ))
3735 bordFace = _borderFaces.insert( faceID ).first;
3737 if ( bordFace == _borderFaces.end() || !intFace ) continue;
3739 // get all links of mesh faces on internal geom face sharing nodes on edge <s>
3740 set< SMESH_OrientedLink > links; //!< links of faces on internal geom face
3741 list<const SMDS_MeshElement*> suspectFaces[2]; //!< mesh faces on border geom faces
3742 int nbSuspectFaces = 0;
3743 SMESHDS_SubMesh* intFaceSM = meshDS->MeshElements( intFace );
3744 if ( !intFaceSM || intFaceSM->NbElements() == 0 ) continue;
3745 SMESH_subMeshIteratorPtr smIt = mesh.GetSubMesh( s )->getDependsOnIterator(true,true);
3746 while ( smIt->more() )
3748 SMESHDS_SubMesh* sm = smIt->next()->GetSubMeshDS();
3749 if ( !sm ) continue;
3750 SMDS_NodeIteratorPtr nIt = sm->GetNodes();
3751 while ( nIt->more() )
3753 const SMDS_MeshNode* nOnEdge = nIt->next();
3754 SMDS_ElemIteratorPtr fIt = nOnEdge->GetInverseElementIterator(SMDSAbs_Face);
3755 while ( fIt->more() )
3757 const SMDS_MeshElement* f = fIt->next();
3758 const int nbNodes = f->NbCornerNodes();
3759 if ( intFaceSM->Contains( f ))
3761 for ( int i = 0; i < nbNodes; ++i )
3762 links.insert( SMESH_OrientedLink( f->GetNode(i), f->GetNode((i+1)%nbNodes)));
3767 for ( int i = 0; i < nbNodes; ++i )
3768 nbDblNodes += isInternalShape( f->GetNode(i)->getshapeId() );
3770 suspectFaces[ nbDblNodes < 2 ].push_back( f );
3776 // suspectFaces[0] having link with same orientation as mesh faces on
3777 // the internal geom face are <borderElems>. suspectFaces[1] have
3778 // only one node on edge <s>, we decide on them later (at the 2nd loop)
3779 // by links of <borderElems> found at the 1st and 2nd loops
3780 set< SMESH_OrientedLink > borderLinks;
3781 for ( int isPostponed = 0; isPostponed < 2; ++isPostponed )
3783 list<const SMDS_MeshElement*>::iterator fIt = suspectFaces[isPostponed].begin();
3784 for ( int nbF = 0; fIt != suspectFaces[isPostponed].end(); ++fIt, ++nbF )
3786 const SMDS_MeshElement* f = *fIt;
3787 bool isBorder = false, linkFound = false, borderLinkFound = false;
3788 list< SMESH_OrientedLink > faceLinks;
3789 int nbNodes = f->NbCornerNodes();
3790 for ( int i = 0; i < nbNodes; ++i )
3792 SMESH_OrientedLink link( f->GetNode(i), f->GetNode((i+1)%nbNodes));
3793 faceLinks.push_back( link );
3796 set< SMESH_OrientedLink >::iterator foundLink = links.find( link );
3797 if ( foundLink != links.end() )
3800 isBorder = ( foundLink->_reversed == link._reversed );
3801 if ( !isBorder && !isPostponed ) break;
3802 faceLinks.pop_back();
3804 else if ( isPostponed && !borderLinkFound )
3806 foundLink = borderLinks.find( link );
3807 if ( foundLink != borderLinks.end() )
3809 borderLinkFound = true;
3810 isBorder = ( foundLink->_reversed != link._reversed );
3817 borderElems.insert( f );
3818 borderLinks.insert( faceLinks.begin(), faceLinks.end() );
3820 else if ( !linkFound && !borderLinkFound )
3822 suspectFaces[1].push_back( f );
3823 if ( nbF > 2 * nbSuspectFaces )
3824 break; // dead loop protection
3831 //================================================================================
3833 * \brief put internal shapes in maps and fill in submeshes to precompute
3835 //================================================================================
3837 void NETGENPlugin_Internals::getInternalEdges( TopTools_IndexedMapOfShape& fmap,
3838 TopTools_IndexedMapOfShape& emap,
3839 TopTools_IndexedMapOfShape& vmap,
3840 list< SMESH_subMesh* > smToPrecompute[])
3842 if ( !hasInternalEdges() ) return;
3843 map<int,int>::const_iterator ev_face = _e2face.begin();
3844 for ( ; ev_face != _e2face.end(); ++ev_face )
3846 const TopoDS_Shape& ev = _mesh.GetMeshDS()->IndexToShape( ev_face->first );
3847 const TopoDS_Shape& face = _mesh.GetMeshDS()->IndexToShape( ev_face->second );
3849 ( ev.ShapeType() == TopAbs_EDGE ? emap : vmap ).Add( ev );
3851 //cout<<"INTERNAL EDGE or VERTEX "<<ev_face->first<<" on face "<<ev_face->second<<endl;
3853 smToPrecompute[ MeshDim_1D ].push_back( _mesh.GetSubMeshContaining( ev_face->first ));
3857 //================================================================================
3859 * \brief return shapes and submeshes to be meshed and already meshed boundary submeshes
3861 //================================================================================
3863 void NETGENPlugin_Internals::getInternalFaces( TopTools_IndexedMapOfShape& fmap,
3864 TopTools_IndexedMapOfShape& emap,
3865 list< SMESH_subMesh* >& intFaceSM,
3866 list< SMESH_subMesh* >& boundarySM)
3868 if ( !hasInternalFaces() ) return;
3870 // <fmap> and <emap> are for not yet meshed shapes
3871 // <intFaceSM> is for submeshes of faces
3872 // <boundarySM> is for meshed edges and vertices
3877 set<int> shapeIDs ( _intShapes );
3878 if ( !_borderFaces.empty() )
3879 shapeIDs.insert( _borderFaces.begin(), _borderFaces.end() );
3881 set<int>::const_iterator intS = shapeIDs.begin();
3882 for ( ; intS != shapeIDs.end(); ++intS )
3884 SMESH_subMesh* sm = _mesh.GetSubMeshContaining( *intS );
3886 if ( sm->GetSubShape().ShapeType() != TopAbs_FACE ) continue;
3888 intFaceSM.push_back( sm );
3890 // add submeshes of not computed internal faces
3891 if ( !sm->IsEmpty() ) continue;
3893 SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(true,true);
3894 while ( smIt->more() )
3897 const TopoDS_Shape& s = sm->GetSubShape();
3899 if ( sm->IsEmpty() )
3902 switch ( s.ShapeType() ) {
3903 case TopAbs_FACE: fmap.Add ( s ); break;
3904 case TopAbs_EDGE: emap.Add ( s ); break;
3910 if ( s.ShapeType() != TopAbs_FACE )
3911 boundarySM.push_back( sm );
3917 //================================================================================
3919 * \brief Return true if given shape is to be precomputed in order to be correctly
3920 * added to netgen mesh
3922 //================================================================================
3924 bool NETGENPlugin_Internals::isShapeToPrecompute(const TopoDS_Shape& s)
3926 int shapeID = _mesh.GetMeshDS()->ShapeToIndex( s );
3927 switch ( s.ShapeType() ) {
3928 case TopAbs_FACE : break; //return isInternalShape( shapeID ) || isBorderFace( shapeID );
3929 case TopAbs_EDGE : return isInternalEdge( shapeID );
3930 case TopAbs_VERTEX: break;
3936 //================================================================================
3938 * \brief Return SMESH
3940 //================================================================================
3942 SMESH_Mesh& NETGENPlugin_Internals::getMesh() const
3944 return const_cast<SMESH_Mesh&>( _mesh );
3947 //================================================================================
3949 * \brief Access to a counter of NETGENPlugin_NetgenLibWrapper instances
3951 //================================================================================
3953 int& NETGENPlugin_NetgenLibWrapper::instanceCounter()
3955 static int theCouner = 0;
3959 //================================================================================
3961 * \brief Initialize netgen library
3963 //================================================================================
3965 NETGENPlugin_NetgenLibWrapper::NETGENPlugin_NetgenLibWrapper()
3967 if ( instanceCounter() == 0 )
3970 ++instanceCounter();
3972 _isComputeOk = false;
3976 if ( !getenv( "KEEP_NETGEN_OUTPUT" ))
3978 // redirect all netgen output (mycout,myerr,cout) to _outputFileName
3979 _outputFileName = getOutputFileName();
3980 _ngcout = netgen::mycout;
3981 _ngcerr = netgen::myerr;
3982 netgen::mycout = new ofstream ( _outputFileName.c_str() );
3983 netgen::myerr = netgen::mycout;
3984 _coutBuffer = std::cout.rdbuf();
3986 cout << "NOTE: netgen output is redirected to file " << _outputFileName << endl;
3988 std::cout.rdbuf( netgen::mycout->rdbuf() );
3992 _ngMesh = Ng_NewMesh();
3995 //================================================================================
3997 * \brief Finish using netgen library
3999 //================================================================================
4001 NETGENPlugin_NetgenLibWrapper::~NETGENPlugin_NetgenLibWrapper()
4003 --instanceCounter();
4005 Ng_DeleteMesh( _ngMesh );
4009 std::cout.rdbuf( _coutBuffer );
4016 //================================================================================
4018 * \brief Set netgen mesh to delete at destruction
4020 //================================================================================
4022 void NETGENPlugin_NetgenLibWrapper::setMesh( Ng_Mesh* mesh )
4025 Ng_DeleteMesh( _ngMesh );
4029 //================================================================================
4031 * \brief Return a unique file name
4033 //================================================================================
4035 std::string NETGENPlugin_NetgenLibWrapper::getOutputFileName()
4037 std::string aTmpDir = SALOMEDS_Tool::GetTmpDir();
4039 TCollection_AsciiString aGenericName = (char*)aTmpDir.c_str();
4040 aGenericName += "NETGEN_";
4042 aGenericName += getpid();
4044 aGenericName += _getpid();
4046 aGenericName += "_";
4047 aGenericName += Abs((Standard_Integer)(long) aGenericName.ToCString());
4048 aGenericName += ".out";
4050 return aGenericName.ToCString();
4053 //================================================================================
4055 * \brief Remove "test.out" and "problemfaces" files in current directory
4057 //================================================================================
4059 void NETGENPlugin_NetgenLibWrapper::RemoveTmpFiles()
4061 bool rm = SMESH_File("test.out").remove() ;
4063 if ( rm && netgen::testout && instanceCounter() == 0 )
4065 delete netgen::testout;
4066 netgen::testout = 0;
4069 SMESH_File("problemfaces").remove();
4070 SMESH_File("occmesh.rep").remove();
4073 //================================================================================
4075 * \brief Remove file with netgen output
4077 //================================================================================
4079 void NETGENPlugin_NetgenLibWrapper::removeOutputFile()
4081 if ( !_outputFileName.empty() )
4085 delete netgen::mycout;
4086 netgen::mycout = _ngcout;
4087 netgen::myerr = _ngcerr;
4090 string tmpDir = SALOMEDS_Tool::GetDirFromPath ( _outputFileName );
4091 string aFileName = SALOMEDS_Tool::GetNameFromPath( _outputFileName ) + ".out";
4092 SALOMEDS::ListOfFileNames_var aFiles = new SALOMEDS::ListOfFileNames;
4094 aFiles[0] = aFileName.c_str();
4096 SALOMEDS_Tool::RemoveTemporaryFiles( tmpDir.c_str(), aFiles.in(), true );