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 //=============================================================================
130 //=============================================================================
132 NETGENPlugin_Mesher::NETGENPlugin_Mesher (SMESH_Mesh* mesh,
133 const TopoDS_Shape& aShape,
139 _fineness(NETGENPlugin_Hypothesis::GetDefaultFineness()),
140 _isViscousLayers2D(false),
149 SetDefaultParameters();
150 ShapesWithLocalSize.Clear();
151 VertexId2LocalSize.clear();
152 EdgeId2LocalSize.clear();
153 FaceId2LocalSize.clear();
154 SolidId2LocalSize.clear();
157 //================================================================================
161 //================================================================================
163 NETGENPlugin_Mesher::~NETGENPlugin_Mesher()
171 //================================================================================
173 * Set pointer to NETGENPlugin_Mesher* field of the holder, that will be
174 * nullified at destruction of this
176 //================================================================================
178 void NETGENPlugin_Mesher::SetSelfPointer( NETGENPlugin_Mesher ** ptr )
189 //================================================================================
191 * \brief Initialize global NETGEN parameters with default values
193 //================================================================================
195 void NETGENPlugin_Mesher::SetDefaultParameters()
197 netgen::MeshingParameters& mparams = netgen::mparam;
198 // maximal mesh edge size
199 mparams.maxh = 0;//NETGENPlugin_Hypothesis::GetDefaultMaxSize();
201 // minimal number of segments per edge
202 mparams.segmentsperedge = NETGENPlugin_Hypothesis::GetDefaultNbSegPerEdge();
203 // rate of growth of size between elements
204 mparams.grading = NETGENPlugin_Hypothesis::GetDefaultGrowthRate();
205 // safety factor for curvatures (elements per radius)
206 mparams.curvaturesafety = NETGENPlugin_Hypothesis::GetDefaultNbSegPerRadius();
207 // create elements of second order
208 mparams.secondorder = NETGENPlugin_Hypothesis::GetDefaultSecondOrder();
209 // quad-dominated surface meshing
213 mparams.quad = NETGENPlugin_Hypothesis_2D::GetDefaultQuadAllowed();
214 _fineness = NETGENPlugin_Hypothesis::GetDefaultFineness();
215 mparams.uselocalh = NETGENPlugin_Hypothesis::GetDefaultSurfaceCurvature();
216 netgen::merge_solids = NETGENPlugin_Hypothesis::GetDefaultFuseEdges();
219 //=============================================================================
223 //=============================================================================
225 void SetLocalSize(TopoDS_Shape GeomShape, double LocalSize)
227 if ( GeomShape.IsNull() ) return;
228 TopAbs_ShapeEnum GeomType = GeomShape.ShapeType();
229 if (GeomType == TopAbs_COMPOUND) {
230 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()) {
231 SetLocalSize(it.Value(), LocalSize);
236 if (! ShapesWithLocalSize.Contains(GeomShape))
237 key = ShapesWithLocalSize.Add(GeomShape);
239 key = ShapesWithLocalSize.FindIndex(GeomShape);
240 if (GeomType == TopAbs_VERTEX) {
241 VertexId2LocalSize[key] = LocalSize;
242 } else if (GeomType == TopAbs_EDGE) {
243 EdgeId2LocalSize[key] = LocalSize;
244 } else if (GeomType == TopAbs_FACE) {
245 FaceId2LocalSize[key] = LocalSize;
246 } else if (GeomType == TopAbs_SOLID) {
247 SolidId2LocalSize[key] = LocalSize;
251 //=============================================================================
253 * Pass parameters to NETGEN
255 //=============================================================================
256 void NETGENPlugin_Mesher::SetParameters(const NETGENPlugin_Hypothesis* hyp)
260 netgen::MeshingParameters& mparams = netgen::mparam;
261 // Initialize global NETGEN parameters:
262 // maximal mesh segment size
263 mparams.maxh = hyp->GetMaxSize();
264 // maximal mesh element linear size
265 mparams.minh = hyp->GetMinSize();
266 // minimal number of segments per edge
267 mparams.segmentsperedge = hyp->GetNbSegPerEdge();
268 // rate of growth of size between elements
269 mparams.grading = hyp->GetGrowthRate();
270 // safety factor for curvatures (elements per radius)
271 mparams.curvaturesafety = hyp->GetNbSegPerRadius();
272 // create elements of second order
273 mparams.secondorder = hyp->GetSecondOrder() ? 1 : 0;
274 // quad-dominated surface meshing
275 mparams.quad = hyp->GetQuadAllowed() ? 1 : 0;
276 _optimize = hyp->GetOptimize();
277 _fineness = hyp->GetFineness();
278 mparams.uselocalh = hyp->GetSurfaceCurvature();
279 netgen::merge_solids = hyp->GetFuseEdges();
282 SMESH_Gen_i* smeshGen_i = SMESH_Gen_i::GetSMESHGen();
283 CORBA::Object_var anObject = smeshGen_i->GetNS()->Resolve("/myStudyManager");
284 SALOMEDS::StudyManager_var aStudyMgr = SALOMEDS::StudyManager::_narrow(anObject);
285 SALOMEDS::Study_var myStudy = aStudyMgr->GetStudyByID(hyp->GetStudyId());
287 const NETGENPlugin_Hypothesis::TLocalSize localSizes = hyp->GetLocalSizesAndEntries();
288 NETGENPlugin_Hypothesis::TLocalSize::const_iterator it = localSizes.begin();
289 for ( ; it != localSizes.end() ; it++)
291 std::string entry = (*it).first;
292 double val = (*it).second;
294 GEOM::GEOM_Object_var aGeomObj;
295 SALOMEDS::SObject_var aSObj = myStudy->FindObjectID( entry.c_str() );
296 if ( !aSObj->_is_nil() ) {
297 CORBA::Object_var obj = aSObj->GetObject();
298 aGeomObj = GEOM::GEOM_Object::_narrow(obj);
301 TopoDS_Shape S = smeshGen_i->GeomObjectToShape( aGeomObj.in() );
302 SetLocalSize(S, val);
307 //=============================================================================
309 * Pass simple parameters to NETGEN
311 //=============================================================================
313 void NETGENPlugin_Mesher::SetParameters(const NETGENPlugin_SimpleHypothesis_2D* hyp)
317 SetDefaultParameters();
320 //=============================================================================
322 * Link - a pair of integer numbers
324 //=============================================================================
328 Link(int _n1, int _n2) : n1(_n1), n2(_n2) {}
329 Link() : n1(0), n2(0) {}
330 bool Contains( int n ) const { return n == n1 || n == n2; }
331 bool IsConnected( const Link& other ) const
333 return (( Contains( other.n1 ) || Contains( other.n2 )) && ( this != &other ));
337 int HashCode(const Link& aLink, int aLimit)
339 return HashCode(aLink.n1 + aLink.n2, aLimit);
342 Standard_Boolean IsEqual(const Link& aLink1, const Link& aLink2)
344 return (( aLink1.n1 == aLink2.n1 && aLink1.n2 == aLink2.n2 ) ||
345 ( aLink1.n1 == aLink2.n2 && aLink1.n2 == aLink2.n1 ));
350 //================================================================================
352 * \brief return id of netgen point corresponding to SMDS node
354 //================================================================================
355 typedef map< const SMDS_MeshNode*, int > TNode2IdMap;
357 int ngNodeId( const SMDS_MeshNode* node,
358 netgen::Mesh& ngMesh,
359 TNode2IdMap& nodeNgIdMap)
361 int newNgId = ngMesh.GetNP() + 1;
363 TNode2IdMap::iterator node_id = nodeNgIdMap.insert( make_pair( node, newNgId )).first;
365 if ( node_id->second == newNgId)
367 #if defined(DUMP_SEGMENTS) || defined(DUMP_TRIANGLES)
368 cout << "Ng " << newNgId << " - " << node;
370 netgen::MeshPoint p( netgen::Point<3> (node->X(), node->Y(), node->Z()) );
371 ngMesh.AddPoint( p );
373 return node_id->second;
376 //================================================================================
378 * \brief Return computed EDGEs connected to the given one
380 //================================================================================
382 list< TopoDS_Edge > getConnectedEdges( const TopoDS_Edge& edge,
383 const TopoDS_Face& face,
384 const set< SMESH_subMesh* > & computedSM,
385 const SMESH_MesherHelper& helper,
386 map< SMESH_subMesh*, set< int > >& addedEdgeSM2Faces)
389 list< TopoDS_Edge > edges;
390 list< int > nbEdgesInWire;
391 /*int nbWires =*/ SMESH_Block::GetOrderedEdges( face, edges, nbEdgesInWire);
393 // find <edge> within <edges>
394 list< TopoDS_Edge >::iterator eItFwd = edges.begin();
395 for ( ; eItFwd != edges.end(); ++eItFwd )
396 if ( edge.IsSame( *eItFwd ))
398 if ( eItFwd == edges.end()) return list< TopoDS_Edge>();
400 if ( eItFwd->Orientation() >= TopAbs_INTERNAL )
402 // connected INTERNAL edges returned from GetOrderedEdges() are wrongly oriented
403 // so treat each INTERNAL edge separately
404 TopoDS_Edge e = *eItFwd;
406 edges.push_back( e );
410 // get all computed EDGEs connected to <edge>
412 list< TopoDS_Edge >::iterator eItBack = eItFwd, ePrev;
413 TopoDS_Vertex vCommon;
414 TopTools_MapOfShape eAdded; // map used not to add a seam edge twice to <edges>
417 // put edges before <edge> to <edges> back
418 while ( edges.begin() != eItFwd )
419 edges.splice( edges.end(), edges, edges.begin() );
423 while ( ++eItFwd != edges.end() )
425 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( *eItFwd );
427 bool connected = TopExp::CommonVertex( *ePrev, *eItFwd, vCommon );
428 bool computed = sm->IsMeshComputed();
429 bool added = addedEdgeSM2Faces[ sm ].count( helper.GetSubShapeID() );
430 bool doubled = !eAdded.Add( *eItFwd );
431 bool orientOK = (( ePrev ->Orientation() < TopAbs_INTERNAL ) ==
432 ( eItFwd->Orientation() < TopAbs_INTERNAL ) );
433 if ( !connected || !computed || !orientOK || added || doubled )
435 // stop advancement; move edges from tail to head
436 while ( edges.back() != *ePrev )
437 edges.splice( edges.begin(), edges, --edges.end() );
443 while ( eItBack != edges.begin() )
447 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( *eItBack );
449 bool connected = TopExp::CommonVertex( *ePrev, *eItBack, vCommon );
450 bool computed = sm->IsMeshComputed();
451 bool added = addedEdgeSM2Faces[ sm ].count( helper.GetSubShapeID() );
452 bool doubled = !eAdded.Add( *eItBack );
453 bool orientOK = (( ePrev ->Orientation() < TopAbs_INTERNAL ) ==
454 ( eItBack->Orientation() < TopAbs_INTERNAL ) );
455 if ( !connected || !computed || !orientOK || added || doubled)
458 edges.erase( edges.begin(), ePrev );
462 if ( edges.front() != edges.back() )
464 // assure that the 1st vertex is meshed
465 TopoDS_Edge eLast = edges.back();
466 while ( !SMESH_Algo::VertexNode( SMESH_MesherHelper::IthVertex( 0, edges.front()), helper.GetMeshDS())
468 edges.front() != eLast )
469 edges.splice( edges.end(), edges, edges.begin() );
474 //================================================================================
476 * \brief Make triangulation of a shape precise enough
478 //================================================================================
480 void updateTriangulation( const TopoDS_Shape& shape )
482 // static set< Poly_Triangulation* > updated;
484 // TopLoc_Location loc;
485 // TopExp_Explorer fExp( shape, TopAbs_FACE );
486 // for ( ; fExp.More(); fExp.Next() )
488 // Handle(Poly_Triangulation) triangulation =
489 // BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
490 // if ( triangulation.IsNull() ||
491 // updated.insert( triangulation.operator->() ).second )
493 // BRepTools::Clean (shape);
496 BRepMesh_IncrementalMesh e(shape, 0.01, true);
498 catch (Standard_Failure)
501 // updated.erase( triangulation.operator->() );
502 // triangulation = BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
503 // updated.insert( triangulation.operator->() );
507 //================================================================================
509 * \brief Returns a medium node either existing in SMESH of created by NETGEN
510 * \param [in] corner1 - corner node 1
511 * \param [in] corner2 - corner node 2
512 * \param [in] defaultMedium - the node created by NETGEN
513 * \param [in] helper - holder of medium nodes existing in SMESH
514 * \return const SMDS_MeshNode* - the result node
516 //================================================================================
518 const SMDS_MeshNode* mediumNode( const SMDS_MeshNode* corner1,
519 const SMDS_MeshNode* corner2,
520 const SMDS_MeshNode* defaultMedium,
521 const SMESH_MesherHelper* helper)
525 TLinkNodeMap::const_iterator l2n =
526 helper->GetTLinkNodeMap().find( SMESH_TLink( corner1, corner2 ));
527 if ( l2n != helper->GetTLinkNodeMap().end() )
528 defaultMedium = l2n->second;
530 return defaultMedium;
533 //================================================================================
535 * \brief Assure that mesh on given shapes is quadratic
537 //================================================================================
539 // void makeQuadratic( const TopTools_IndexedMapOfShape& shapes,
540 // SMESH_Mesh* mesh )
542 // for ( int i = 1; i <= shapes.Extent(); ++i )
544 // SMESHDS_SubMesh* smDS = mesh->GetMeshDS()->MeshElements( shapes(i) );
545 // if ( !smDS ) continue;
546 // SMDS_ElemIteratorPtr elemIt = smDS->GetElements();
547 // if ( !elemIt->more() ) continue;
548 // const SMDS_MeshElement* e = elemIt->next();
549 // if ( !e || e->IsQuadratic() )
552 // TIDSortedElemSet elems;
553 // elems.insert( e );
554 // while ( elemIt->more() )
555 // elems.insert( elems.end(), elemIt->next() );
557 // SMESH_MeshEditor( mesh ).ConvertToQuadratic( /*3d=*/false, elems, /*biQuad=*/false );
563 //================================================================================
565 * \brief Initialize netgen::OCCGeometry with OCCT shape
567 //================================================================================
569 void NETGENPlugin_Mesher::PrepareOCCgeometry(netgen::OCCGeometry& occgeo,
570 const TopoDS_Shape& shape,
572 list< SMESH_subMesh* > * meshedSM,
573 NETGENPlugin_Internals* intern)
575 updateTriangulation( shape );
578 BRepBndLib::Add (shape, bb);
579 double x1,y1,z1,x2,y2,z2;
580 bb.Get (x1,y1,z1,x2,y2,z2);
581 netgen::Point<3> p1 = netgen::Point<3> (x1,y1,z1);
582 netgen::Point<3> p2 = netgen::Point<3> (x2,y2,z2);
583 occgeo.boundingbox = netgen::Box<3> (p1,p2);
585 occgeo.shape = shape;
588 // fill maps of shapes of occgeo with not yet meshed subshapes
590 // get root submeshes
591 list< SMESH_subMesh* > rootSM;
592 const int shapeID = mesh.GetMeshDS()->ShapeToIndex( shape );
593 if ( shapeID > 0 ) { // SMESH_subMesh with ID 0 may exist, don't use it!
594 rootSM.push_back( mesh.GetSubMesh( shape ));
597 for ( TopoDS_Iterator it( shape ); it.More(); it.Next() )
598 rootSM.push_back( mesh.GetSubMesh( it.Value() ));
603 // add subshapes of empty submeshes
604 list< SMESH_subMesh* >::iterator rootIt = rootSM.begin(), rootEnd = rootSM.end();
605 for ( ; rootIt != rootEnd; ++rootIt ) {
606 SMESH_subMesh * root = *rootIt;
607 SMESH_subMeshIteratorPtr smIt = root->getDependsOnIterator(/*includeSelf=*/true,
608 /*complexShapeFirst=*/true);
609 // to find a right orientation of subshapes (PAL20462)
610 TopTools_IndexedMapOfShape subShapes;
611 TopExp::MapShapes(root->GetSubShape(), subShapes);
612 while ( smIt->more() )
614 SMESH_subMesh* sm = smIt->next();
615 TopoDS_Shape shape = sm->GetSubShape();
616 totNbFaces += ( shape.ShapeType() == TopAbs_FACE );
617 if ( intern && intern->isShapeToPrecompute( shape ))
619 if ( !meshedSM || sm->IsEmpty() )
621 if ( shape.ShapeType() != TopAbs_VERTEX )
622 shape = subShapes( subShapes.FindIndex( shape ));// shape -> index -> oriented shape
623 if ( shape.Orientation() >= TopAbs_INTERNAL )
624 shape.Orientation( TopAbs_FORWARD ); // isuue 0020676
625 switch ( shape.ShapeType() ) {
626 case TopAbs_FACE : occgeo.fmap.Add( shape ); break;
627 case TopAbs_EDGE : occgeo.emap.Add( shape ); break;
628 case TopAbs_VERTEX: occgeo.vmap.Add( shape ); break;
629 case TopAbs_SOLID :occgeo.somap.Add( shape ); break;
633 // collect submeshes of meshed shapes
636 const int dim = SMESH_Gen::GetShapeDim( shape );
637 meshedSM[ dim ].push_back( sm );
641 occgeo.facemeshstatus.SetSize (totNbFaces);
642 occgeo.facemeshstatus = 0;
643 occgeo.face_maxh_modified.SetSize(totNbFaces);
644 occgeo.face_maxh_modified = 0;
645 occgeo.face_maxh.SetSize(totNbFaces);
646 occgeo.face_maxh = netgen::mparam.maxh;
649 //================================================================================
651 * \brief Return a default min size value suitable for the given geometry.
653 //================================================================================
655 double NETGENPlugin_Mesher::GetDefaultMinSize(const TopoDS_Shape& geom,
656 const double maxSize)
658 updateTriangulation( geom );
662 const int* pi[4] = { &i1, &i2, &i3, &i1 };
665 TopExp_Explorer fExp( geom, TopAbs_FACE );
666 for ( ; fExp.More(); fExp.Next() )
668 Handle(Poly_Triangulation) triangulation =
669 BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
670 if ( triangulation.IsNull() ) continue;
671 const double fTol = BRep_Tool::Tolerance( TopoDS::Face( fExp.Current() ));
672 const TColgp_Array1OfPnt& points = triangulation->Nodes();
673 const Poly_Array1OfTriangle& trias = triangulation->Triangles();
674 for ( int iT = trias.Lower(); iT <= trias.Upper(); ++iT )
676 trias(iT).Get( i1, i2, i3 );
677 for ( int j = 0; j < 3; ++j )
679 double dist2 = points(*pi[j]).SquareDistance( points( *pi[j+1] ));
680 if ( dist2 < minh && fTol*fTol < dist2 )
682 bb.Add( points(*pi[j]));
686 if ( minh > 0.25 * bb.SquareExtent() ) // simple geometry, rough triangulation
688 minh = 1e-3 * sqrt( bb.SquareExtent());
689 //cout << "BND BOX minh = " <<minh << endl;
693 minh = 3 * sqrt( minh ); // triangulation for visualization is rather fine
694 //cout << "TRIANGULATION minh = " <<minh << endl;
696 if ( minh > 0.5 * maxSize )
702 //================================================================================
704 * \brief Restrict size of elements at a given point
706 //================================================================================
708 void NETGENPlugin_Mesher::RestrictLocalSize(netgen::Mesh& ngMesh,
711 const bool overrideMinH)
713 if ( size <= std::numeric_limits<double>::min() )
715 if ( netgen::mparam.minh > size )
719 ngMesh.SetMinimalH( size );
720 netgen::mparam.minh = size;
724 size = netgen::mparam.minh;
727 netgen::Point3d pi(p.X(), p.Y(), p.Z());
728 ngMesh.RestrictLocalH( pi, size );
731 //================================================================================
733 * \brief fill ngMesh with nodes and elements of computed submeshes
735 //================================================================================
737 bool NETGENPlugin_Mesher::FillNgMesh(netgen::OCCGeometry& occgeom,
738 netgen::Mesh& ngMesh,
739 vector<const SMDS_MeshNode*>& nodeVec,
740 const list< SMESH_subMesh* > & meshedSM,
741 SMESH_MesherHelper* quadHelper,
742 SMESH_ProxyMesh::Ptr proxyMesh)
744 TNode2IdMap nodeNgIdMap;
745 for ( size_t i = 1; i < nodeVec.size(); ++i )
746 nodeNgIdMap.insert( make_pair( nodeVec[i], i ));
748 TopTools_MapOfShape visitedShapes;
749 map< SMESH_subMesh*, set< int > > visitedEdgeSM2Faces;
750 set< SMESH_subMesh* > computedSM( meshedSM.begin(), meshedSM.end() );
752 SMESH_MesherHelper helper (*_mesh);
754 int faceNgID = ngMesh.GetNFD();
756 list< SMESH_subMesh* >::const_iterator smIt, smEnd = meshedSM.end();
757 for ( smIt = meshedSM.begin(); smIt != smEnd; ++smIt )
759 SMESH_subMesh* sm = *smIt;
760 if ( !visitedShapes.Add( sm->GetSubShape() ))
763 const SMESHDS_SubMesh * smDS = sm->GetSubMeshDS();
764 if ( !smDS ) continue;
766 switch ( sm->GetSubShape().ShapeType() )
768 case TopAbs_EDGE: { // EDGE
769 // ----------------------
770 TopoDS_Edge geomEdge = TopoDS::Edge( sm->GetSubShape() );
771 if ( geomEdge.Orientation() >= TopAbs_INTERNAL )
772 geomEdge.Orientation( TopAbs_FORWARD ); // issue 0020676
774 // Add ng segments for each not meshed FACE the EDGE bounds
775 PShapeIteratorPtr fIt = helper.GetAncestors( geomEdge, *sm->GetFather(), TopAbs_FACE );
776 while ( const TopoDS_Shape * anc = fIt->next() )
778 faceNgID = occgeom.fmap.FindIndex( *anc );
780 continue; // meshed face
782 int faceSMDSId = helper.GetMeshDS()->ShapeToIndex( *anc );
783 if ( visitedEdgeSM2Faces[ sm ].count( faceSMDSId ))
784 continue; // already treated EDGE
786 TopoDS_Face face = TopoDS::Face( occgeom.fmap( faceNgID ));
787 if ( face.Orientation() >= TopAbs_INTERNAL )
788 face.Orientation( TopAbs_FORWARD ); // issue 0020676
790 // get all meshed EDGEs of the FACE connected to geomEdge (issue 0021140)
791 helper.SetSubShape( face );
792 list< TopoDS_Edge > edges = getConnectedEdges( geomEdge, face, computedSM, helper,
793 visitedEdgeSM2Faces );
795 continue; // wrong ancestor?
797 // find out orientation of <edges> within <face>
798 TopoDS_Edge eNotSeam = edges.front();
799 if ( helper.HasSeam() )
801 list< TopoDS_Edge >::iterator eIt = edges.begin();
802 while ( helper.IsRealSeam( *eIt )) ++eIt;
803 if ( eIt != edges.end() )
806 TopAbs_Orientation fOri = helper.GetSubShapeOri( face, eNotSeam );
807 bool isForwad = ( fOri == eNotSeam.Orientation() || fOri >= TopAbs_INTERNAL );
809 // get all nodes from connected <edges>
810 const bool isQuad = smDS->IsQuadratic();
811 StdMeshers_FaceSide fSide( face, edges, _mesh, isForwad, isQuad );
812 const vector<UVPtStruct>& points = fSide.GetUVPtStruct();
813 if ( points.empty() )
814 return false; // invalid node params?
815 int i, nbSeg = fSide.NbSegments();
817 // remember EDGEs of fSide to treat only once
818 for ( int iE = 0; iE < fSide.NbEdges(); ++iE )
819 visitedEdgeSM2Faces[ helper.GetMesh()->GetSubMesh( fSide.Edge(iE )) ].insert(faceSMDSId);
821 double otherSeamParam = 0;
826 int prevNgId = ngNodeId( points[0].node, ngMesh, nodeNgIdMap );
828 for ( i = 0; i < nbSeg; ++i )
830 const UVPtStruct& p1 = points[ i ];
831 const UVPtStruct& p2 = points[ i+1 ];
833 if ( p1.node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX ) //an EDGE begins
836 if ( helper.IsRealSeam( p1.node->getshapeId() ))
838 TopoDS_Edge e = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
839 isSeam = helper.IsRealSeam( e );
842 otherSeamParam = helper.GetOtherParam( helper.GetPeriodicIndex() & 1 ? p2.u : p2.v );
849 seg[1] = prevNgId = ngNodeId( p2.node, ngMesh, nodeNgIdMap );
850 // node param on curve
851 seg.epgeominfo[ 0 ].dist = p1.param;
852 seg.epgeominfo[ 1 ].dist = p2.param;
854 seg.epgeominfo[ 0 ].u = p1.u;
855 seg.epgeominfo[ 0 ].v = p1.v;
856 seg.epgeominfo[ 1 ].u = p2.u;
857 seg.epgeominfo[ 1 ].v = p2.v;
859 //geomEdge = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
860 //seg.epgeominfo[ 0 ].edgenr = seg.epgeominfo[ 1 ].edgenr = occgeom.emap.FindIndex( geomEdge );
862 //seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
863 seg.si = faceNgID; // = geom.fmap.FindIndex (face);
864 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
865 ngMesh.AddSegment (seg);
867 SMESH_TNodeXYZ np1( p1.node ), np2( p2.node );
868 RestrictLocalSize( ngMesh, 0.5*(np1+np2), (np1-np2).Modulus() );
871 cout << "Segment: " << seg.edgenr << " on SMESH face " << helper.GetMeshDS()->ShapeToIndex( face ) << endl
872 << "\tface index: " << seg.si << endl
873 << "\tp1: " << seg[0] << endl
874 << "\tp2: " << seg[1] << endl
875 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
876 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
877 //<< "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
878 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
879 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl;
880 //<< "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
884 if ( helper.GetPeriodicIndex() && 1 ) {
885 seg.epgeominfo[ 0 ].u = otherSeamParam;
886 seg.epgeominfo[ 1 ].u = otherSeamParam;
887 swap (seg.epgeominfo[0].v, seg.epgeominfo[1].v);
889 seg.epgeominfo[ 0 ].v = otherSeamParam;
890 seg.epgeominfo[ 1 ].v = otherSeamParam;
891 swap (seg.epgeominfo[0].u, seg.epgeominfo[1].u);
893 swap( seg[0], seg[1] );
894 swap( seg.epgeominfo[0].dist, seg.epgeominfo[1].dist );
895 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
896 ngMesh.AddSegment( seg );
898 cout << "Segment: " << seg.edgenr << endl
899 << "\t is SEAM (reverse) of the previous. "
900 << " Other " << (helper.GetPeriodicIndex() && 1 ? "U" : "V")
901 << " = " << otherSeamParam << endl;
904 else if ( fOri == TopAbs_INTERNAL )
906 swap( seg[0], seg[1] );
907 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
908 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
909 ngMesh.AddSegment( seg );
911 cout << "Segment: " << seg.edgenr << endl << "\t is REVERSE of the previous" << endl;
915 } // loop on geomEdge ancestors
917 if ( quadHelper ) // remember medium nodes of sub-meshes
919 SMDS_ElemIteratorPtr edges = smDS->GetElements();
920 while ( edges->more() )
922 const SMDS_MeshElement* e = edges->next();
923 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshEdge*>( e )))
929 } // case TopAbs_EDGE
931 case TopAbs_FACE: { // FACE
932 // ----------------------
933 const TopoDS_Face& geomFace = TopoDS::Face( sm->GetSubShape() );
934 helper.SetSubShape( geomFace );
935 bool isInternalFace = ( geomFace.Orientation() == TopAbs_INTERNAL );
937 // Find solids the geomFace bounds
938 int solidID1 = 0, solidID2 = 0;
939 StdMeshers_QuadToTriaAdaptor* quadAdaptor =
940 dynamic_cast<StdMeshers_QuadToTriaAdaptor*>( proxyMesh.get() );
943 solidID1 = occgeom.somap.FindIndex( quadAdaptor->GetShape() );
947 PShapeIteratorPtr solidIt = helper.GetAncestors( geomFace, *sm->GetFather(), TopAbs_SOLID);
948 while ( const TopoDS_Shape * solid = solidIt->next() )
950 int id = occgeom.somap.FindIndex ( *solid );
951 if ( solidID1 && id != solidID1 ) solidID2 = id;
955 // Add ng face descriptors of meshed faces
957 ngMesh.AddFaceDescriptor( netgen::FaceDescriptor( faceNgID, solidID1, solidID2, 0 ));
959 // if second oreder is required, even already meshed faces must be passed to NETGEN
960 int fID = occgeom.fmap.Add( geomFace );
961 if ( occgeom.facemeshstatus.Size() < fID ) occgeom.facemeshstatus.SetSize( fID );
962 occgeom.facemeshstatus[ fID-1 ] = netgen::FACE_MESHED_OK;
963 while ( fID < faceNgID ) // geomFace is already in occgeom.fmap, add a copy
965 fID = occgeom.fmap.Add( BRepBuilderAPI_Copy( geomFace, /*copyGeom=*/false ));
966 if ( occgeom.facemeshstatus.Size() < fID ) occgeom.facemeshstatus.SetSize( fID );
967 occgeom.facemeshstatus[ fID-1 ] = netgen::FACE_MESHED_OK;
969 // Problem with the second order in a quadrangular mesh remains.
970 // 1) All quadrangles generated by NETGEN are moved to an inexistent face
971 // by FillSMesh() (find "AddFaceDescriptor")
972 // 2) Temporary triangles generated by StdMeshers_QuadToTriaAdaptor
973 // are on faces where quadrangles were.
974 // Due to these 2 points, wrong geom faces are used while conversion to quadratic
975 // of the mentioned above quadrangles and triangles
977 // Orient the face correctly in solidID1 (issue 0020206)
978 bool reverse = false;
980 TopoDS_Shape solid = occgeom.somap( solidID1 );
981 TopAbs_Orientation faceOriInSolid = helper.GetSubShapeOri( solid, geomFace );
982 if ( faceOriInSolid >= 0 )
984 helper.IsReversedSubMesh( TopoDS::Face( geomFace.Oriented( faceOriInSolid )));
987 // Add surface elements
989 netgen::Element2d tri(3);
990 tri.SetIndex( faceNgID );
991 SMESH_TNodeXYZ xyz[3];
993 #ifdef DUMP_TRIANGLES
994 cout << "SMESH face " << helper.GetMeshDS()->ShapeToIndex( geomFace )
995 << " internal="<<isInternalFace << endl;
998 smDS = proxyMesh->GetSubMesh( geomFace );
1000 SMDS_ElemIteratorPtr faces = smDS->GetElements();
1001 while ( faces->more() )
1003 const SMDS_MeshElement* f = faces->next();
1004 if ( f->NbNodes() % 3 != 0 ) // not triangle
1006 PShapeIteratorPtr solidIt=helper.GetAncestors(geomFace,*sm->GetFather(),TopAbs_SOLID);
1007 if ( const TopoDS_Shape * solid = solidIt->next() )
1008 sm = _mesh->GetSubMesh( *solid );
1009 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1010 smError.reset( new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH,"Not triangle sub-mesh"));
1011 smError->myBadElements.push_back( f );
1015 for ( int i = 0; i < 3; ++i )
1017 const SMDS_MeshNode* node = f->GetNode( i ), * inFaceNode=0;
1020 // get node UV on face
1021 int shapeID = node->getshapeId();
1022 if ( helper.IsSeamShape( shapeID ))
1024 if ( helper.IsSeamShape( f->GetNodeWrap( i+1 )->getshapeId() ))
1025 inFaceNode = f->GetNodeWrap( i-1 );
1027 inFaceNode = f->GetNodeWrap( i+1 );
1029 gp_XY uv = helper.GetNodeUV( geomFace, node, inFaceNode );
1031 int ind = reverse ? 3-i : i+1;
1032 tri.GeomInfoPi(ind).u = uv.X();
1033 tri.GeomInfoPi(ind).v = uv.Y();
1034 tri.PNum (ind) = ngNodeId( node, ngMesh, nodeNgIdMap );
1037 // pass a triangle size to NG size-map
1038 double size = ( ( xyz[0] - xyz[1] ).Modulus() +
1039 ( xyz[1] - xyz[2] ).Modulus() +
1040 ( xyz[2] - xyz[0] ).Modulus() ) / 3;
1041 gp_XYZ gc = ( xyz[0] + xyz[1] + xyz[2] ) / 3;
1042 RestrictLocalSize( ngMesh, gc, size, /*overrideMinH=*/false );
1044 ngMesh.AddSurfaceElement (tri);
1045 #ifdef DUMP_TRIANGLES
1046 cout << tri << endl;
1049 if ( isInternalFace )
1051 swap( tri[1], tri[2] );
1052 ngMesh.AddSurfaceElement (tri);
1053 #ifdef DUMP_TRIANGLES
1054 cout << tri << endl;
1059 if ( quadHelper ) // remember medium nodes of sub-meshes
1061 SMDS_ElemIteratorPtr faces = smDS->GetElements();
1062 while ( faces->more() )
1064 const SMDS_MeshElement* f = faces->next();
1065 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshFace*>( f )))
1071 } // case TopAbs_FACE
1073 case TopAbs_VERTEX: { // VERTEX
1074 // --------------------------
1075 // issue 0021405. Add node only if a VERTEX is shared by a not meshed EDGE,
1076 // else netgen removes a free node and nodeVector becomes invalid
1077 PShapeIteratorPtr ansIt = helper.GetAncestors( sm->GetSubShape(),
1081 while ( const TopoDS_Shape* e = ansIt->next() )
1083 SMESH_subMesh* eSub = helper.GetMesh()->GetSubMesh( *e );
1084 if (( toAdd = ( eSub->IsEmpty() && !SMESH_Algo::isDegenerated( TopoDS::Edge( *e )))))
1089 SMDS_NodeIteratorPtr nodeIt = smDS->GetNodes();
1090 if ( nodeIt->more() )
1091 ngNodeId( nodeIt->next(), ngMesh, nodeNgIdMap );
1097 } // loop on submeshes
1100 nodeVec.resize( ngMesh.GetNP() + 1 );
1101 TNode2IdMap::iterator node_NgId, nodeNgIdEnd = nodeNgIdMap.end();
1102 for ( node_NgId = nodeNgIdMap.begin(); node_NgId != nodeNgIdEnd; ++node_NgId)
1103 nodeVec[ node_NgId->second ] = node_NgId->first;
1108 //================================================================================
1110 * \brief Duplicate mesh faces on internal geom faces
1112 //================================================================================
1114 void NETGENPlugin_Mesher::FixIntFaces(const netgen::OCCGeometry& occgeom,
1115 netgen::Mesh& ngMesh,
1116 NETGENPlugin_Internals& internalShapes)
1118 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1120 // find ng indices of internal faces
1122 for ( int ngFaceID = 1; ngFaceID <= occgeom.fmap.Extent(); ++ngFaceID )
1124 int smeshID = meshDS->ShapeToIndex( occgeom.fmap( ngFaceID ));
1125 if ( internalShapes.isInternalShape( smeshID ))
1126 ngFaceIds.insert( ngFaceID );
1128 if ( !ngFaceIds.empty() )
1131 int i, nbFaces = ngMesh.GetNSE();
1132 for ( i = 1; i <= nbFaces; ++i)
1134 netgen::Element2d elem = ngMesh.SurfaceElement(i);
1135 if ( ngFaceIds.count( elem.GetIndex() ))
1137 swap( elem[1], elem[2] );
1138 ngMesh.AddSurfaceElement (elem);
1144 //================================================================================
1146 * \brief Tries to heal the mesh on a FACE. The FACE is supposed to be partially
1147 * meshed due to NETGEN failure
1148 * \param [in] occgeom - geometry
1149 * \param [in,out] ngMesh - the mesh to fix
1150 * \param [inout] faceID - ID of the FACE to fix the mesh on
1151 * \return bool - is mesh is or becomes OK
1153 //================================================================================
1155 bool NETGENPlugin_Mesher::FixFaceMesh(const netgen::OCCGeometry& occgeom,
1156 netgen::Mesh& ngMesh,
1159 // we address a case where the FACE is almost fully meshed except small holes
1160 // of usually triangular shape at FACE boundary (IPAL52861)
1162 // The case appeared to be not simple: holes only look triangular but
1163 // indeed are a self intersecting polygon. A reason of the bug was in coincident
1164 // NG points on a seam edge. But the code below is very nice, leave it for
1169 if ( occgeom.fmap.Extent() < faceID )
1171 //const TopoDS_Face& face = TopoDS::Face( occgeom.fmap( faceID ));
1173 // find free links on the FACE
1174 NCollection_Map<Link> linkMap;
1175 for ( int iF = 1; iF <= ngMesh.GetNSE(); ++iF )
1177 const netgen::Element2d& elem = ngMesh.SurfaceElement(iF);
1178 if ( faceID != elem.GetIndex() )
1180 int n0 = elem[ elem.GetNP() - 1 ];
1181 for ( int i = 0; i < elem.GetNP(); ++i )
1184 Link link( n0, n1 );
1185 if ( !linkMap.Add( link ))
1186 linkMap.Remove( link );
1190 // add/remove boundary links
1191 for ( int iSeg = 1; iSeg <= ngMesh.GetNSeg(); ++iSeg )
1193 const netgen::Segment& seg = ngMesh.LineSegment( iSeg );
1194 if ( seg.si != faceID ) // !edgeIDs.Contains( seg.edgenr ))
1196 Link link( seg[1], seg[0] ); // reverse!!!
1197 if ( !linkMap.Add( link ))
1198 linkMap.Remove( link );
1200 if ( linkMap.IsEmpty() )
1202 if ( linkMap.Extent() < 3 )
1205 // make triangles of the links
1207 netgen::Element2d tri(3);
1208 tri.SetIndex ( faceID );
1210 NCollection_Map<Link>::Iterator linkIt( linkMap );
1211 Link link1 = linkIt.Value();
1212 // look for a link connected to link1
1213 NCollection_Map<Link>::Iterator linkIt2 = linkIt;
1214 for ( linkIt2.Next(); linkIt2.More(); linkIt2.Next() )
1216 const Link& link2 = linkIt2.Value();
1217 if ( link2.IsConnected( link1 ))
1219 // look for a link connected to both link1 and link2
1220 NCollection_Map<Link>::Iterator linkIt3 = linkIt2;
1221 for ( linkIt3.Next(); linkIt3.More(); linkIt3.Next() )
1223 const Link& link3 = linkIt3.Value();
1224 if ( link3.IsConnected( link1 ) &&
1225 link3.IsConnected( link2 ) )
1230 tri[2] = ( link2.Contains( link1.n1 ) ? link2.n1 : link3.n1 );
1231 if ( tri[0] == tri[2] || tri[1] == tri[2] )
1233 ngMesh.AddSurfaceElement( tri );
1235 // prepare for the next tria search
1236 if ( linkMap.Extent() == 3 )
1238 linkMap.Remove( link3 );
1239 linkMap.Remove( link2 );
1241 linkMap.Remove( link1 );
1242 link1 = linkIt.Value();
1255 //================================================================================
1256 // define gp_XY_Subtracted pointer to function calling gp_XY::Subtracted(gp_XY)
1257 gp_XY_FunPtr(Subtracted);
1258 //gp_XY_FunPtr(Added);
1260 //================================================================================
1262 * \brief Evaluate distance between two 2d points along the surface
1264 //================================================================================
1266 double evalDist( const gp_XY& uv1,
1268 const Handle(Geom_Surface)& surf,
1269 const int stopHandler=-1)
1271 if ( stopHandler > 0 ) // continue recursion
1273 gp_XY mid = SMESH_MesherHelper::GetMiddleUV( surf, uv1, uv2 );
1274 return evalDist( uv1,mid, surf, stopHandler-1 ) + evalDist( mid,uv2, surf, stopHandler-1 );
1276 double dist3D = surf->Value( uv1.X(), uv1.Y() ).Distance( surf->Value( uv2.X(), uv2.Y() ));
1277 if ( stopHandler == 0 ) // stop recursion
1280 // start recursion if necessary
1281 double dist2D = SMESH_MesherHelper::ApplyIn2D(surf, uv1, uv2, gp_XY_Subtracted, 0).Modulus();
1282 if ( fabs( dist3D - dist2D ) < dist2D * 1e-10 )
1283 return dist3D; // equal parametrization of a planar surface
1285 return evalDist( uv1, uv2, surf, 3 ); // start recursion
1288 //================================================================================
1290 * \brief Data of vertex internal in geom face
1292 //================================================================================
1296 gp_XY uv; //!< UV in face parametric space
1297 int ngId; //!< ng id of corrsponding node
1298 gp_XY uvClose; //!< UV of closest boundary node
1299 int ngIdClose; //!< ng id of closest boundary node
1302 //================================================================================
1304 * \brief Data of vertex internal in solid
1306 //================================================================================
1310 int ngId; //!< ng id of corresponding node
1311 int ngIdClose; //!< ng id of closest 2d mesh element
1312 int ngIdCloseN; //!< ng id of closest node of the closest 2d mesh element
1315 inline double dist2(const netgen::MeshPoint& p1, const netgen::MeshPoint& p2)
1317 return gp_Pnt( NGPOINT_COORDS(p1)).SquareDistance( gp_Pnt( NGPOINT_COORDS(p2)));
1321 //================================================================================
1323 * \brief Make netgen take internal vertices in faces into account by adding
1324 * segments including internal vertices
1326 * This function works in supposition that 1D mesh is already computed in ngMesh
1328 //================================================================================
1330 void NETGENPlugin_Mesher::AddIntVerticesInFaces(const netgen::OCCGeometry& occgeom,
1331 netgen::Mesh& ngMesh,
1332 vector<const SMDS_MeshNode*>& nodeVec,
1333 NETGENPlugin_Internals& internalShapes)
1335 if ((int) nodeVec.size() < ngMesh.GetNP() )
1336 nodeVec.resize( ngMesh.GetNP(), 0 );
1338 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1339 SMESH_MesherHelper helper( internalShapes.getMesh() );
1341 const map<int,list<int> >& face2Vert = internalShapes.getFacesWithVertices();
1342 map<int,list<int> >::const_iterator f2v = face2Vert.begin();
1343 for ( ; f2v != face2Vert.end(); ++f2v )
1345 const TopoDS_Face& face = TopoDS::Face( meshDS->IndexToShape( f2v->first ));
1346 if ( face.IsNull() ) continue;
1347 int faceNgID = occgeom.fmap.FindIndex (face);
1348 if ( faceNgID < 0 ) continue;
1350 TopLoc_Location loc;
1351 Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc);
1353 helper.SetSubShape( face );
1354 helper.SetElementsOnShape( true );
1356 // Get data of internal vertices and add them to ngMesh
1358 multimap< double, TIntVData > dist2VData; // sort vertices by distance from boundary nodes
1360 int i, nbSegInit = ngMesh.GetNSeg();
1362 // boundary characteristics
1363 double totSegLen2D = 0;
1366 const list<int>& iVertices = f2v->second;
1367 list<int>::const_iterator iv = iVertices.begin();
1368 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1371 // get node on vertex
1372 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1373 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1376 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1377 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1378 nV = SMESH_Algo::VertexNode( V, meshDS );
1379 if ( !nV ) continue;
1382 netgen::MeshPoint mp( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1383 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1384 vData.ngId = ngMesh.GetNP();
1385 nodeVec.push_back( nV );
1389 vData.uv = helper.GetNodeUV( face, nV, 0, &uvOK );
1390 if ( !uvOK ) helper.CheckNodeUV( face, nV, vData.uv, BRep_Tool::Tolerance(V),/*force=*/1);
1392 // loop on all segments of the face to find the node closest to vertex and to count
1393 // average segment 2d length
1394 double closeDist2 = numeric_limits<double>::max(), dist2;
1396 for (i = 1; i <= ngMesh.GetNSeg(); ++i)
1398 netgen::Segment & seg = ngMesh.LineSegment(i);
1399 if ( seg.si != faceNgID ) continue;
1401 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1403 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1404 if ( ngIdLast == seg[ iEnd ] ) continue;
1405 dist2 = helper.ApplyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1406 if ( dist2 < closeDist2 )
1407 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1408 ngIdLast = seg[ iEnd ];
1412 totSegLen2D += helper.ApplyIn2D(surf, uv[0], uv[1], gp_XY_Subtracted, false).Modulus();
1416 dist2VData.insert( make_pair( closeDist2, vData ));
1419 if ( totNbSeg == 0 ) break;
1420 double avgSegLen2d = totSegLen2D / totNbSeg;
1422 // Loop on vertices to add segments
1424 multimap< double, TIntVData >::iterator dist_vData = dist2VData.begin();
1425 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1427 double closeDist2 = dist_vData->first, dist2;
1428 TIntVData & vData = dist_vData->second;
1430 // try to find more close node among segments added for internal vertices
1431 for (i = nbSegInit+1; i <= ngMesh.GetNSeg(); ++i)
1433 netgen::Segment & seg = ngMesh.LineSegment(i);
1434 if ( seg.si != faceNgID ) continue;
1436 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1438 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1439 dist2 = helper.ApplyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1440 if ( dist2 < closeDist2 )
1441 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1444 // decide whether to use the closest node as the second end of segment or to
1445 // create a new point
1446 int segEnd1 = vData.ngId;
1447 int segEnd2 = vData.ngIdClose; // to use closest node
1448 gp_XY uvV = vData.uv, uvP = vData.uvClose;
1449 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1450 double nodeDist2D = sqrt( closeDist2 );
1451 double nodeDist3D = evalDist( vData.uv, vData.uvClose, surf );
1452 bool avgLenOK = ( avgSegLen2d < 0.75 * nodeDist2D );
1453 bool hintLenOK = ( segLenHint < 0.75 * nodeDist3D );
1454 //cout << "uvV " << uvV.X() <<","<<uvV.Y() << " ";
1455 if ( hintLenOK || avgLenOK )
1457 // create a point between the closest node and V
1460 double r = min( 0.5, ( hintLenOK ? segLenHint/nodeDist3D : avgSegLen2d/nodeDist2D ));
1461 // direction from V to closet node in 2D
1462 gp_Dir2d v2n( helper.ApplyIn2D(surf, uvP, uvV, gp_XY_Subtracted, false ));
1464 uvP = vData.uv + r * nodeDist2D * v2n.XY();
1465 gp_Pnt P = surf->Value( uvP.X(), uvP.Y() ).Transformed( loc );
1467 netgen::MeshPoint mp( netgen::Point<3> (P.X(), P.Y(), P.Z()));
1468 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1469 segEnd2 = ngMesh.GetNP();
1470 //cout << "Middle " << r << " uv " << uvP.X() << "," << uvP.Y() << "( " << ngMesh.Point(segEnd2).X()<<","<<ngMesh.Point(segEnd2).Y()<<","<<ngMesh.Point(segEnd2).Z()<<" )"<< endl;
1471 SMDS_MeshNode * nP = helper.AddNode(P.X(), P.Y(), P.Z());
1472 nodeVec.push_back( nP );
1474 //else cout << "at Node " << " uv " << uvP.X() << "," << uvP.Y() << endl;
1477 netgen::Segment seg;
1479 if ( segEnd1 > segEnd2 ) swap( segEnd1, segEnd2 ), swap( uvV, uvP );
1480 seg[0] = segEnd1; // ng node id
1481 seg[1] = segEnd2; // ng node id
1482 seg.edgenr = ngMesh.GetNSeg() + 1;// segment id
1485 seg.epgeominfo[ 0 ].dist = 0; // param on curve
1486 seg.epgeominfo[ 0 ].u = uvV.X();
1487 seg.epgeominfo[ 0 ].v = uvV.Y();
1488 seg.epgeominfo[ 1 ].dist = 1; // param on curve
1489 seg.epgeominfo[ 1 ].u = uvP.X();
1490 seg.epgeominfo[ 1 ].v = uvP.Y();
1492 // seg.epgeominfo[ 0 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1493 // seg.epgeominfo[ 1 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1495 ngMesh.AddSegment (seg);
1497 // add reverse segment
1498 swap( seg[0], seg[1] );
1499 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1500 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1501 ngMesh.AddSegment (seg);
1507 //================================================================================
1509 * \brief Make netgen take internal vertices in solids into account by adding
1510 * faces including internal vertices
1512 * This function works in supposition that 2D mesh is already computed in ngMesh
1514 //================================================================================
1516 void NETGENPlugin_Mesher::AddIntVerticesInSolids(const netgen::OCCGeometry& occgeom,
1517 netgen::Mesh& ngMesh,
1518 vector<const SMDS_MeshNode*>& nodeVec,
1519 NETGENPlugin_Internals& internalShapes)
1521 #ifdef DUMP_TRIANGLES_SCRIPT
1522 // create a python script making a mesh containing triangles added for internal vertices
1523 ofstream py(DUMP_TRIANGLES_SCRIPT);
1524 py << "import SMESH"<< endl
1525 << "from salome.smesh import smeshBuilder"<<endl
1526 << "smesh = smeshBuilder.New(salome.myStudy)"<<endl
1527 << "m = smesh.Mesh(name='triangles')" << endl;
1529 if ((int) nodeVec.size() < ngMesh.GetNP() )
1530 nodeVec.resize( ngMesh.GetNP(), 0 );
1532 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1533 SMESH_MesherHelper helper( internalShapes.getMesh() );
1535 const map<int,list<int> >& so2Vert = internalShapes.getSolidsWithVertices();
1536 map<int,list<int> >::const_iterator s2v = so2Vert.begin();
1537 for ( ; s2v != so2Vert.end(); ++s2v )
1539 const TopoDS_Shape& solid = meshDS->IndexToShape( s2v->first );
1540 if ( solid.IsNull() ) continue;
1541 int solidNgID = occgeom.somap.FindIndex (solid);
1542 if ( solidNgID < 0 && !occgeom.somap.IsEmpty() ) continue;
1544 helper.SetSubShape( solid );
1545 helper.SetElementsOnShape( true );
1547 // find ng indices of faces within the solid
1549 for (TopExp_Explorer fExp(solid, TopAbs_FACE); fExp.More(); fExp.Next() )
1550 ngFaceIds.insert( occgeom.fmap.FindIndex( fExp.Current() ));
1551 if ( ngFaceIds.size() == 1 && *ngFaceIds.begin() == 0 )
1552 ngFaceIds.insert( 1 );
1554 // Get data of internal vertices and add them to ngMesh
1556 multimap< double, TIntVSoData > dist2VData; // sort vertices by distance from ng faces
1558 int i, nbFaceInit = ngMesh.GetNSE();
1560 // boundary characteristics
1561 double totSegLen = 0;
1564 const list<int>& iVertices = s2v->second;
1565 list<int>::const_iterator iv = iVertices.begin();
1566 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1569 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1571 // get node on vertex
1572 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1575 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1576 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1577 nV = SMESH_Algo::VertexNode( V, meshDS );
1578 if ( !nV ) continue;
1581 netgen::MeshPoint mpV( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1582 ngMesh.AddPoint ( mpV, 1, netgen::FIXEDPOINT );
1583 vData.ngId = ngMesh.GetNP();
1584 nodeVec.push_back( nV );
1586 // loop on all 2d elements to find the one closest to vertex and to count
1587 // average segment length
1588 double closeDist2 = numeric_limits<double>::max(), avgDist2;
1589 for (i = 1; i <= ngMesh.GetNSE(); ++i)
1591 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1592 if ( !ngFaceIds.count( elem.GetIndex() )) continue;
1594 multimap< double, int> dist2nID; // sort nodes of element by distance from V
1595 for ( int j = 0; j < elem.GetNP(); ++j)
1597 netgen::MeshPoint mp = ngMesh.Point( elem[j] );
1598 double d2 = dist2( mpV, mp );
1599 dist2nID.insert( make_pair( d2, elem[j] ));
1600 avgDist2 += d2 / elem.GetNP();
1602 totNbSeg++, totSegLen+= sqrt( dist2( mp, ngMesh.Point( elem[(j+1)%elem.GetNP()])));
1604 double dist = dist2nID.begin()->first; //avgDist2;
1605 if ( dist < closeDist2 )
1606 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= dist;
1608 dist2VData.insert( make_pair( closeDist2, vData ));
1611 if ( totNbSeg == 0 ) break;
1612 double avgSegLen = totSegLen / totNbSeg;
1614 // Loop on vertices to add triangles
1616 multimap< double, TIntVSoData >::iterator dist_vData = dist2VData.begin();
1617 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1619 double closeDist2 = dist_vData->first;
1620 TIntVSoData & vData = dist_vData->second;
1622 const netgen::MeshPoint& mpV = ngMesh.Point( vData.ngId );
1624 // try to find more close face among ones added for internal vertices
1625 for (i = nbFaceInit+1; i <= ngMesh.GetNSE(); ++i)
1627 double avgDist2 = 0;
1628 multimap< double, int> dist2nID;
1629 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1630 for ( int j = 0; j < elem.GetNP(); ++j)
1632 double d = dist2( mpV, ngMesh.Point( elem[j] ));
1633 dist2nID.insert( make_pair( d, elem[j] ));
1634 avgDist2 += d / elem.GetNP();
1635 if ( avgDist2 < closeDist2 )
1636 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= avgDist2;
1639 // sort nodes of the closest face by angle with vector from V to the closest node
1640 const double tol = numeric_limits<double>::min();
1641 map< double, int > angle2ID;
1642 const netgen::Element2d& closeFace = ngMesh.SurfaceElement( vData.ngIdClose );
1643 netgen::MeshPoint mp[2];
1644 mp[0] = ngMesh.Point( vData.ngIdCloseN );
1645 gp_XYZ p1( NGPOINT_COORDS( mp[0] ));
1646 gp_XYZ pV( NGPOINT_COORDS( mpV ));
1647 gp_Vec v2p1( pV, p1 );
1648 double distN1 = v2p1.Magnitude();
1649 if ( distN1 <= tol ) continue;
1651 for ( int j = 0; j < closeFace.GetNP(); ++j)
1653 mp[1] = ngMesh.Point( closeFace[j] );
1654 gp_Vec v2p( pV, gp_Pnt( NGPOINT_COORDS( mp[1] )) );
1655 angle2ID.insert( make_pair( v2p1.Angle( v2p ), closeFace[j]));
1657 // get node with angle of 60 degrees or greater
1658 map< double, int >::iterator angle_id = angle2ID.lower_bound( 60. * M_PI / 180. );
1659 if ( angle_id == angle2ID.end() ) angle_id = --angle2ID.end();
1660 const double minAngle = 30. * M_PI / 180.;
1661 const double angle = angle_id->first;
1662 bool angleOK = ( angle > minAngle );
1664 // find points to create a triangle
1665 netgen::Element2d tri(3);
1667 tri[0] = vData.ngId;
1668 tri[1] = vData.ngIdCloseN; // to use the closest nodes
1669 tri[2] = angle_id->second; // to use the node with best angle
1671 // decide whether to use the closest node and the node with best angle or to create new ones
1672 for ( int isBestAngleN = 0; isBestAngleN < 2; ++isBestAngleN )
1674 bool createNew = !angleOK; //, distOK = true;
1676 int triInd = isBestAngleN ? 2 : 1;
1677 mp[isBestAngleN] = ngMesh.Point( tri[triInd] );
1682 double distN2 = sqrt( dist2( mpV, mp[isBestAngleN]));
1683 createNew = ( fabs( distN2 - distN1 ) > 0.25 * distN1 );
1685 else if ( angle < tol )
1687 v2p1.SetX( v2p1.X() + 1e-3 );
1693 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1694 bool avgLenOK = ( avgSegLen < 0.75 * distN1 );
1695 bool hintLenOK = ( segLenHint < 0.75 * distN1 );
1696 createNew = (createNew || avgLenOK || hintLenOK );
1697 // we create a new node not closer than 0.5 to the closest face
1698 // in order not to clash with other close face
1699 double r = min( 0.5, ( hintLenOK ? segLenHint : avgSegLen ) / distN1 );
1700 distFromV = r * distN1;
1704 // create a new point, between the node and the vertex if angleOK
1705 gp_XYZ p( NGPOINT_COORDS( mp[isBestAngleN] ));
1706 gp_Vec v2p( pV, p ); v2p.Normalize();
1707 if ( isBestAngleN && !angleOK )
1708 p = p1 + gp_Dir( v2p.XYZ() - v2p1.XYZ()).XYZ() * distN1 * 0.95;
1710 p = pV + v2p.XYZ() * distFromV;
1712 if ( !isBestAngleN ) p1 = p, distN1 = distFromV;
1714 mp[isBestAngleN].SetPoint( netgen::Point<3> (p.X(), p.Y(), p.Z()));
1715 ngMesh.AddPoint ( mp[isBestAngleN], 1, netgen::SURFACEPOINT );
1716 tri[triInd] = ngMesh.GetNP();
1717 nodeVec.push_back( helper.AddNode( p.X(), p.Y(), p.Z()) );
1720 ngMesh.AddSurfaceElement (tri);
1721 swap( tri[1], tri[2] );
1722 ngMesh.AddSurfaceElement (tri);
1724 #ifdef DUMP_TRIANGLES_SCRIPT
1725 py << "n1 = m.AddNode( "<< mpV(0)<<", "<< mpV(1)<<", "<< mpV(2)<<") "<< endl
1726 << "n2 = m.AddNode( "<< mp[0](0)<<", "<< mp[0](1)<<", "<< mp[0](2)<<") "<< endl
1727 << "n3 = m.AddNode( "<< mp[1](0)<<", "<< mp[1](1)<<", "<< mp[1](2)<<" )" << endl
1728 << "m.AddFace([n1,n2,n3])" << endl;
1730 } // loop on internal vertices of a solid
1732 } // loop on solids with internal vertices
1735 //================================================================================
1737 * \brief Fill netgen mesh with segments of a FACE
1738 * \param ngMesh - netgen mesh
1739 * \param geom - container of OCCT geometry to mesh
1740 * \param wires - data of nodes on FACE boundary
1741 * \param helper - mesher helper holding the FACE
1742 * \param nodeVec - vector of nodes in which node index == netgen ID
1743 * \retval SMESH_ComputeErrorPtr - error description
1745 //================================================================================
1747 SMESH_ComputeErrorPtr
1748 NETGENPlugin_Mesher::AddSegmentsToMesh(netgen::Mesh& ngMesh,
1749 netgen::OCCGeometry& geom,
1750 const TSideVector& wires,
1751 SMESH_MesherHelper& helper,
1752 vector< const SMDS_MeshNode* > & nodeVec,
1753 const bool overrideMinH)
1755 // ----------------------------
1756 // Check wires and count nodes
1757 // ----------------------------
1759 for ( size_t iW = 0; iW < wires.size(); ++iW )
1761 StdMeshers_FaceSidePtr wire = wires[ iW ];
1762 if ( wire->MissVertexNode() )
1764 // Commented for issue 0020960. It worked for the case, let's wait for case where it doesn't.
1765 // It seems that there is no reason for this limitation
1767 // (new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH, "Missing nodes on vertices"));
1769 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1770 if ((int) uvPtVec.size() != wire->NbPoints() )
1771 return SMESH_ComputeError::New(COMPERR_BAD_INPUT_MESH,
1772 SMESH_Comment("Unexpected nb of points on wire ") << iW
1773 << ": " << uvPtVec.size()<<" != "<<wire->NbPoints());
1774 nbNodes += wire->NbPoints();
1776 nodeVec.reserve( nodeVec.size() + nbNodes + 1 );
1777 if ( nodeVec.empty() )
1778 nodeVec.push_back( 0 );
1780 // -----------------
1782 // -----------------
1784 const bool wasNgMeshEmpty = ( ngMesh.GetNP() < 1 ); /* true => this method is called by
1785 NETGENPlugin_NETGEN_2D_ONLY */
1787 // map for nodes on vertices since they can be shared between wires
1788 // ( issue 0020676, face_int_box.brep) and nodes built by NETGEN
1789 map<const SMDS_MeshNode*, int > node2ngID;
1790 if ( !wasNgMeshEmpty ) // fill node2ngID with nodes built by NETGEN
1792 set< int > subIDs; // ids of sub-shapes of the FACE
1793 for ( size_t iW = 0; iW < wires.size(); ++iW )
1795 StdMeshers_FaceSidePtr wire = wires[ iW ];
1796 for ( int iE = 0, nbE = wire->NbEdges(); iE < nbE; ++iE )
1798 subIDs.insert( wire->EdgeID( iE ));
1799 subIDs.insert( helper.GetMeshDS()->ShapeToIndex( wire->FirstVertex( iE )));
1802 for ( size_t ngID = 1; ngID < nodeVec.size(); ++ngID )
1803 if ( subIDs.count( nodeVec[ngID]->getshapeId() ))
1804 node2ngID.insert( make_pair( nodeVec[ngID], ngID ));
1807 const int solidID = 0, faceID = geom.fmap.FindIndex( helper.GetSubShape() );
1808 if ( ngMesh.GetNFD() < 1 )
1809 ngMesh.AddFaceDescriptor( netgen::FaceDescriptor( faceID, solidID, solidID, 0 ));
1811 for ( size_t iW = 0; iW < wires.size(); ++iW )
1813 StdMeshers_FaceSidePtr wire = wires[ iW ];
1814 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1815 const int nbSegments = wire->NbPoints() - 1;
1817 // assure the 1st node to be in node2ngID, which is needed to correctly
1818 // "close chain of segments" (see below) in case if the 1st node is not
1819 // onVertex because it is on a Viscous layer
1820 node2ngID.insert( make_pair( uvPtVec[ 0 ].node, ngMesh.GetNP() + 1 ));
1822 // compute length of every segment
1823 vector<double> segLen( nbSegments );
1824 for ( int i = 0; i < nbSegments; ++i )
1825 segLen[i] = SMESH_TNodeXYZ( uvPtVec[ i ].node ).Distance( uvPtVec[ i+1 ].node );
1827 int edgeID = 1, posID = -2;
1828 bool isInternalWire = false;
1829 double vertexNormPar = 0;
1830 const int prevNbNGSeg = ngMesh.GetNSeg();
1831 for ( int i = 0; i < nbSegments; ++i ) // loop on segments
1833 // Add the first point of a segment
1835 const SMDS_MeshNode * n = uvPtVec[ i ].node;
1836 const int posShapeID = n->getshapeId();
1837 bool onVertex = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX );
1838 bool onEdge = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_EDGE );
1840 // skip nodes on degenerated edges
1841 if ( helper.IsDegenShape( posShapeID ) &&
1842 helper.IsDegenShape( uvPtVec[ i+1 ].node->getshapeId() ))
1845 int ngID1 = ngMesh.GetNP() + 1, ngID2 = ngID1+1;
1846 if ( onVertex || ( !wasNgMeshEmpty && onEdge ) || helper.IsRealSeam( posShapeID ))
1847 ngID1 = node2ngID.insert( make_pair( n, ngID1 )).first->second;
1848 if ( ngID1 > ngMesh.GetNP() )
1850 netgen::MeshPoint mp( netgen::Point<3> (n->X(), n->Y(), n->Z()) );
1851 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1852 nodeVec.push_back( n );
1854 else // n is in ngMesh already, and ngID2 in prev segment is wrong
1856 ngID2 = ngMesh.GetNP() + 1;
1857 if ( i > 0 ) // prev segment belongs to same wire
1859 netgen::Segment& prevSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
1866 netgen::Segment seg;
1868 seg[0] = ngID1; // ng node id
1869 seg[1] = ngID2; // ng node id
1870 seg.edgenr = ngMesh.GetNSeg() + 1; // ng segment id
1871 seg.si = faceID; // = geom.fmap.FindIndex (face);
1873 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1875 const UVPtStruct& pnt = uvPtVec[ i + iEnd ];
1877 seg.epgeominfo[ iEnd ].dist = pnt.param; // param on curve
1878 seg.epgeominfo[ iEnd ].u = pnt.u;
1879 seg.epgeominfo[ iEnd ].v = pnt.v;
1881 // find out edge id and node parameter on edge
1882 onVertex = ( pnt.normParam + 1e-10 > vertexNormPar );
1883 if ( onVertex || posShapeID != posID )
1886 double normParam = pnt.normParam;
1888 normParam = 0.5 * ( uvPtVec[ i ].normParam + uvPtVec[ i+1 ].normParam );
1889 int edgeIndexInWire = wire->EdgeIndex( normParam );
1890 vertexNormPar = wire->LastParameter( edgeIndexInWire );
1891 const TopoDS_Edge& edge = wire->Edge( edgeIndexInWire );
1892 edgeID = geom.emap.FindIndex( edge );
1894 isInternalWire = ( edge.Orientation() == TopAbs_INTERNAL );
1895 // if ( onVertex ) // param on curve is different on each of two edges
1896 // seg.epgeominfo[ iEnd ].dist = helper.GetNodeU( edge, pnt.node );
1898 seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
1901 ngMesh.AddSegment (seg);
1903 // restrict size of elements near the segment
1904 SMESH_TNodeXYZ np1( n ), np2( uvPtVec[ i+1 ].node );
1905 // get an average size of adjacent segments to avoid sharp change of
1906 // element size (regression on issue 0020452, note 0010898)
1907 int iPrev = SMESH_MesherHelper::WrapIndex( i-1, nbSegments );
1908 int iNext = SMESH_MesherHelper::WrapIndex( i+1, nbSegments );
1909 double sumH = segLen[ iPrev ] + segLen[ i ] + segLen[ iNext ];
1910 int nbSeg = ( int( segLen[ iPrev ] > sumH / 100.) +
1911 int( segLen[ i ] > sumH / 100.) +
1912 int( segLen[ iNext ] > sumH / 100.));
1914 RestrictLocalSize( ngMesh, 0.5*(np1+np2), sumH / nbSeg, overrideMinH );
1916 if ( isInternalWire )
1918 swap (seg[0], seg[1]);
1919 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1920 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1921 ngMesh.AddSegment (seg);
1923 } // loop on segments on a wire
1925 // close chain of segments
1926 if ( nbSegments > 0 )
1928 netgen::Segment& lastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() - int( isInternalWire ));
1929 const SMDS_MeshNode * lastNode = uvPtVec.back().node;
1930 lastSeg[1] = node2ngID.insert( make_pair( lastNode, lastSeg[1] )).first->second;
1931 if ( lastSeg[1] > ngMesh.GetNP() )
1933 netgen::MeshPoint mp( netgen::Point<3> (lastNode->X(), lastNode->Y(), lastNode->Z()) );
1934 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1935 nodeVec.push_back( lastNode );
1937 if ( isInternalWire )
1939 netgen::Segment& realLastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
1940 realLastSeg[0] = lastSeg[1];
1944 #ifdef DUMP_SEGMENTS
1945 cout << "BEGIN WIRE " << iW << endl;
1946 for ( int i = prevNbNGSeg+1; i <= ngMesh.GetNSeg(); ++i )
1948 netgen::Segment& seg = ngMesh.LineSegment( i );
1950 netgen::Segment& prevSeg = ngMesh.LineSegment( i-1 );
1951 if ( seg[0] == prevSeg[1] && seg[1] == prevSeg[0] )
1953 cout << "Segment: " << seg.edgenr << endl << "\tis REVRESE of the previous one" << endl;
1957 cout << "Segment: " << seg.edgenr << endl
1958 << "\tp1: " << seg[0] << " n" << nodeVec[ seg[0]]->GetID() << endl
1959 << "\tp2: " << seg[1] << " n" << nodeVec[ seg[1]]->GetID() << endl
1960 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
1961 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
1962 << "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
1963 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
1964 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl
1965 << "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
1967 cout << "--END WIRE " << iW << endl;
1969 SMESH_Comment __not_unused_variable( prevNbNGSeg );
1972 } // loop on WIREs of a FACE
1974 // add a segment instead of an internal vertex
1975 if ( wasNgMeshEmpty )
1977 NETGENPlugin_Internals intShapes( *helper.GetMesh(), helper.GetSubShape(), /*is3D=*/false );
1978 AddIntVerticesInFaces( geom, ngMesh, nodeVec, intShapes );
1980 ngMesh.CalcSurfacesOfNode();
1985 //================================================================================
1987 * \brief Fill SMESH mesh according to contents of netgen mesh
1988 * \param occgeo - container of OCCT geometry to mesh
1989 * \param ngMesh - netgen mesh
1990 * \param initState - bn of entities in netgen mesh before computing
1991 * \param sMesh - SMESH mesh to fill in
1992 * \param nodeVec - vector of nodes in which node index == netgen ID
1993 * \param comment - returns problem description
1994 * \param quadHelper - holder of medium nodes of sub-meshes
1995 * \retval int - error
1997 //================================================================================
1999 int NETGENPlugin_Mesher::FillSMesh(const netgen::OCCGeometry& occgeo,
2000 netgen::Mesh& ngMesh,
2001 const NETGENPlugin_ngMeshInfo& initState,
2003 std::vector<const SMDS_MeshNode*>& nodeVec,
2004 SMESH_Comment& comment,
2005 SMESH_MesherHelper* quadHelper)
2007 int nbNod = ngMesh.GetNP();
2008 int nbSeg = ngMesh.GetNSeg();
2009 int nbFac = ngMesh.GetNSE();
2010 int nbVol = ngMesh.GetNE();
2012 SMESHDS_Mesh* meshDS = sMesh.GetMeshDS();
2014 // quadHelper is used for either
2015 // 1) making quadratic elements when a lower dimention mesh is loaded
2016 // to SMESH before convertion to quadratic by NETGEN
2017 // 2) sewing of quadratic elements with quadratic elements of sub-meshes
2018 if ( quadHelper && !quadHelper->GetIsQuadratic() && quadHelper->GetTLinkNodeMap().empty() )
2021 // -------------------------------------
2022 // Create and insert nodes into nodeVec
2023 // -------------------------------------
2025 nodeVec.resize( nbNod + 1 );
2026 int i, nbInitNod = initState._nbNodes;
2027 for (i = nbInitNod+1; i <= nbNod; ++i )
2029 const netgen::MeshPoint& ngPoint = ngMesh.Point(i);
2030 SMDS_MeshNode* node = NULL;
2031 TopoDS_Vertex aVert;
2032 // First, netgen creates nodes on vertices in occgeo.vmap,
2033 // so node index corresponds to vertex index
2034 // but (issue 0020776) netgen does not create nodes with equal coordinates
2035 if ( i-nbInitNod <= occgeo.vmap.Extent() )
2037 gp_Pnt p ( NGPOINT_COORDS(ngPoint) );
2038 for (int iV = i-nbInitNod; aVert.IsNull() && iV <= occgeo.vmap.Extent(); ++iV)
2040 aVert = TopoDS::Vertex( occgeo.vmap( iV ));
2041 gp_Pnt pV = BRep_Tool::Pnt( aVert );
2042 if ( p.SquareDistance( pV ) > 1e-20 )
2045 node = const_cast<SMDS_MeshNode*>( SMESH_Algo::VertexNode( aVert, meshDS ));
2048 if (!node) // node not found on vertex
2050 node = meshDS->AddNode( NGPOINT_COORDS( ngPoint ));
2051 if (!aVert.IsNull())
2052 meshDS->SetNodeOnVertex(node, aVert);
2057 // -------------------------------------------
2058 // Create mesh segments along geometric edges
2059 // -------------------------------------------
2061 int nbInitSeg = initState._nbSegments;
2062 for (i = nbInitSeg+1; i <= nbSeg; ++i )
2064 const netgen::Segment& seg = ngMesh.LineSegment(i);
2066 int pinds[3] = { seg.pnums[0], seg.pnums[1], seg.pnums[2] };
2069 for (int j=0; j < 3; ++j)
2071 int pind = pinds[j];
2072 if (pind <= 0 || !nodeVec_ACCESS(pind))
2080 int aGeomEdgeInd = seg.epgeominfo[j].edgenr;
2081 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
2082 aEdge = TopoDS::Edge(occgeo.emap(aGeomEdgeInd));
2084 param = seg.epgeominfo[j].dist;
2087 else // middle point
2089 param = param2 * 0.5;
2091 if (!aEdge.IsNull() && nodeVec_ACCESS(pind)->getshapeId() < 1)
2093 meshDS->SetNodeOnEdge(nodeVec_ACCESS(pind), aEdge, param);
2098 SMDS_MeshEdge* edge = 0;
2099 if (nbp == 2) // second order ?
2101 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1])))
2103 if ( quadHelper ) // final mesh must be quadratic
2104 edge = quadHelper->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
2106 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
2110 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
2111 nodeVec_ACCESS(pinds[2])))
2113 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
2114 nodeVec_ACCESS(pinds[2]));
2118 if ( comment.empty() ) comment << "Cannot create a mesh edge";
2119 MESSAGE("Cannot create a mesh edge");
2120 nbSeg = nbFac = nbVol = 0;
2123 if ( !aEdge.IsNull() && edge->getshapeId() < 1 )
2124 meshDS->SetMeshElementOnShape(edge, aEdge);
2126 else if ( comment.empty() )
2128 comment << "Invalid netgen segment #" << i;
2132 // ----------------------------------------
2133 // Create mesh faces along geometric faces
2134 // ----------------------------------------
2136 int nbInitFac = initState._nbFaces;
2137 int quadFaceID = ngMesh.GetNFD() + 1;
2138 if ( nbInitFac < nbFac )
2139 // add a faces descriptor to exclude qudrangle elements generated by NETGEN
2140 // from computation of 3D mesh
2141 ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(quadFaceID, /*solid1=*/0, /*solid2=*/0, 0));
2143 vector<const SMDS_MeshNode*> nodes;
2144 for (i = nbInitFac+1; i <= nbFac; ++i )
2146 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
2147 const int aGeomFaceInd = elem.GetIndex();
2149 if (aGeomFaceInd > 0 && aGeomFaceInd <= occgeo.fmap.Extent())
2150 aFace = TopoDS::Face(occgeo.fmap(aGeomFaceInd));
2152 for ( int j = 1; j <= elem.GetNP(); ++j )
2154 int pind = elem.PNum(j);
2155 if ( pind < 1 || pind >= (int) nodeVec.size() )
2157 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind))
2159 nodes.push_back( node );
2160 if (!aFace.IsNull() && node->getshapeId() < 1)
2162 const netgen::PointGeomInfo& pgi = elem.GeomInfoPi(j);
2163 meshDS->SetNodeOnFace(node, aFace, pgi.u, pgi.v);
2167 if ((int) nodes.size() != elem.GetNP() )
2169 if ( comment.empty() )
2170 comment << "Invalid netgen 2d element #" << i;
2171 continue; // bad node ids
2173 SMDS_MeshFace* face = NULL;
2174 switch (elem.GetType())
2177 if ( quadHelper ) // final mesh must be quadratic
2178 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2]);
2180 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2]);
2183 if ( quadHelper ) // final mesh must be quadratic
2184 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2186 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2187 // exclude qudrangle elements from computation of 3D mesh
2188 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2191 nodes[5] = mediumNode( nodes[0],nodes[1],nodes[5], quadHelper );
2192 nodes[3] = mediumNode( nodes[1],nodes[2],nodes[3], quadHelper );
2193 nodes[4] = mediumNode( nodes[2],nodes[0],nodes[4], quadHelper );
2194 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[5],nodes[3],nodes[4]);
2197 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2198 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2199 nodes[5] = mediumNode( nodes[2],nodes[3],nodes[5], quadHelper );
2200 nodes[6] = mediumNode( nodes[3],nodes[0],nodes[6], quadHelper );
2201 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3],
2202 nodes[4],nodes[7],nodes[5],nodes[6]);
2203 // exclude qudrangle elements from computation of 3D mesh
2204 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2207 MESSAGE("NETGEN created a face of unexpected type, ignoring");
2212 if ( comment.empty() ) comment << "Cannot create a mesh face";
2213 MESSAGE("Cannot create a mesh face");
2214 nbSeg = nbFac = nbVol = 0;
2217 if ( !aFace.IsNull() )
2218 meshDS->SetMeshElementOnShape( face, aFace );
2221 // ------------------
2222 // Create tetrahedra
2223 // ------------------
2225 for ( i = 1; i <= nbVol; ++i )
2227 const netgen::Element& elem = ngMesh.VolumeElement(i);
2228 int aSolidInd = elem.GetIndex();
2229 TopoDS_Solid aSolid;
2230 if ( aSolidInd > 0 && aSolidInd <= occgeo.somap.Extent() )
2231 aSolid = TopoDS::Solid(occgeo.somap(aSolidInd));
2233 for ( int j = 1; j <= elem.GetNP(); ++j )
2235 int pind = elem.PNum(j);
2236 if ( pind < 1 || pind >= (int)nodeVec.size() )
2238 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind) )
2240 nodes.push_back(node);
2241 if ( !aSolid.IsNull() && node->getshapeId() < 1 )
2242 meshDS->SetNodeInVolume(node, aSolid);
2245 if ((int) nodes.size() != elem.GetNP() )
2247 if ( comment.empty() )
2248 comment << "Invalid netgen 3d element #" << i;
2251 SMDS_MeshVolume* vol = NULL;
2252 switch ( elem.GetType() )
2255 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3]);
2258 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2259 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2260 nodes[5] = mediumNode( nodes[2],nodes[0],nodes[5], quadHelper );
2261 nodes[6] = mediumNode( nodes[0],nodes[3],nodes[6], quadHelper );
2262 nodes[8] = mediumNode( nodes[1],nodes[3],nodes[8], quadHelper );
2263 nodes[9] = mediumNode( nodes[2],nodes[3],nodes[9], quadHelper );
2264 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3],
2265 nodes[4],nodes[7],nodes[5],nodes[6],nodes[8],nodes[9]);
2268 MESSAGE("NETGEN created a volume of unexpected type, ignoring");
2273 if ( comment.empty() ) comment << "Cannot create a mesh volume";
2274 MESSAGE("Cannot create a mesh volume");
2275 nbSeg = nbFac = nbVol = 0;
2278 if (!aSolid.IsNull())
2279 meshDS->SetMeshElementOnShape(vol, aSolid);
2281 return comment.empty() ? 0 : 1;
2286 //================================================================================
2288 * \brief Restrict size of elements on the given edge
2290 //================================================================================
2292 void setLocalSize(const TopoDS_Edge& edge,
2296 if ( size <= std::numeric_limits<double>::min() )
2298 Standard_Real u1, u2;
2299 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, u1, u2);
2300 if ( curve.IsNull() )
2302 TopoDS_Iterator vIt( edge );
2303 if ( !vIt.More() ) return;
2304 gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex( vIt.Value() ));
2305 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
2309 const int nb = (int)( 1.5 * SMESH_Algo::EdgeLength( edge ) / size );
2310 Standard_Real delta = (u2-u1)/nb;
2311 for(int i=0; i<nb; i++)
2313 Standard_Real u = u1 + delta*i;
2314 gp_Pnt p = curve->Value(u);
2315 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
2316 netgen::Point3d pi(p.X(), p.Y(), p.Z());
2317 double resultSize = mesh.GetH(pi);
2318 if ( resultSize - size > 0.1*size )
2319 // netgen does restriction iff oldH/newH > 1.2 (localh.cpp:136)
2320 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), resultSize/1.201 );
2325 //================================================================================
2327 * \brief Convert error into text
2329 //================================================================================
2331 std::string text(int err)
2336 SMESH_Comment("Error in netgen::OCCGenerateMesh() at ") << netgen::multithread.task;
2339 //================================================================================
2341 * \brief Convert exception into text
2343 //================================================================================
2345 std::string text(Standard_Failure& ex)
2347 SMESH_Comment str("Exception in netgen::OCCGenerateMesh()");
2348 str << " at " << netgen::multithread.task
2349 << ": " << ex.DynamicType()->Name();
2350 if ( ex.GetMessageString() && strlen( ex.GetMessageString() ))
2351 str << ": " << ex.GetMessageString();
2354 //================================================================================
2356 * \brief Convert exception into text
2358 //================================================================================
2360 std::string text(netgen::NgException& ex)
2362 SMESH_Comment str("NgException");
2363 if ( strlen( netgen::multithread.task ) > 0 )
2364 str << " at " << netgen::multithread.task;
2365 str << ": " << ex.What();
2369 //================================================================================
2371 * \brief Looks for triangles lying on a SOLID
2373 //================================================================================
2375 bool hasBadElemOnSolid( const list<const SMDS_MeshElement*>& elems,
2376 SMESH_subMesh* solidSM )
2378 TopTools_IndexedMapOfShape solidSubs;
2379 TopExp::MapShapes( solidSM->GetSubShape(), solidSubs );
2380 SMESHDS_Mesh* mesh = solidSM->GetFather()->GetMeshDS();
2382 list<const SMDS_MeshElement*>::const_iterator e = elems.begin();
2383 for ( ; e != elems.end(); ++e )
2385 const SMDS_MeshElement* elem = *e;
2386 // if ( elem->GetType() != SMDSAbs_Face ) -- 23047
2388 int nbNodesOnSolid = 0, nbNodes = elem->NbNodes();
2389 SMDS_NodeIteratorPtr nIt = elem->nodeIterator();
2390 while ( nIt->more() )
2392 const SMDS_MeshNode* n = nIt->next();
2393 const TopoDS_Shape& s = mesh->IndexToShape( n->getshapeId() );
2394 nbNodesOnSolid += ( !s.IsNull() && solidSubs.Contains( s ));
2395 if ( nbNodesOnSolid > 2 ||
2396 nbNodesOnSolid == nbNodes)
2403 const double edgeMeshingTime = 0.001;
2404 const double faceMeshingTime = 0.019;
2405 const double edgeFaceMeshingTime = edgeMeshingTime + faceMeshingTime;
2406 const double faceOptimizTime = 0.06;
2407 const double voluMeshingTime = 0.15;
2408 const double volOptimizeTime = 0.77;
2411 //=============================================================================
2413 * Here we are going to use the NETGEN mesher
2415 //=============================================================================
2417 bool NETGENPlugin_Mesher::Compute()
2419 NETGENPlugin_NetgenLibWrapper ngLib;
2421 netgen::MeshingParameters& mparams = netgen::mparam;
2422 MESSAGE("Compute with:\n"
2423 " max size = " << mparams.maxh << "\n"
2424 " segments per edge = " << mparams.segmentsperedge);
2426 " growth rate = " << mparams.grading << "\n"
2427 " elements per radius = " << mparams.curvaturesafety << "\n"
2428 " second order = " << mparams.secondorder << "\n"
2429 " quad allowed = " << mparams.quad << "\n"
2430 " surface curvature = " << mparams.uselocalh << "\n"
2431 " fuse edges = " << netgen::merge_solids);
2433 SMESH_ComputeErrorPtr error = SMESH_ComputeError::New();
2434 SMESH_MesherHelper quadHelper( *_mesh );
2435 quadHelper.SetIsQuadratic( mparams.secondorder );
2437 static string debugFile = "/tmp/ngMesh.py"; /* to call toPython( _ngMesh, debugFile )
2438 while debugging netgen */
2439 // -------------------------
2440 // Prepare OCC geometry
2441 // -------------------------
2443 netgen::OCCGeometry occgeo;
2444 list< SMESH_subMesh* > meshedSM[3]; // for 0-2 dimensions
2445 NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
2446 PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
2449 _totalTime = edgeFaceMeshingTime;
2451 _totalTime += faceOptimizTime;
2453 _totalTime += voluMeshingTime + ( _optimize ? volOptimizeTime : 0 );
2454 double doneTime = 0;
2457 _curShapeIndex = -1;
2459 // -------------------------
2460 // Generate the mesh
2461 // -------------------------
2464 NETGENPlugin_ngMeshInfo initState; // it remembers size of ng mesh equal to size of Smesh
2466 SMESH_Comment comment;
2469 // vector of nodes in which node index == netgen ID
2470 vector< const SMDS_MeshNode* > nodeVec;
2478 // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
2479 mparams.uselocalh = false;
2480 mparams.grading = 0.8; // not limitited size growth
2482 if ( _simpleHyp->GetNumberOfSegments() )
2484 mparams.maxh = occgeo.boundingbox.Diam();
2487 mparams.maxh = _simpleHyp->GetLocalLength();
2490 if ( mparams.maxh == 0.0 )
2491 mparams.maxh = occgeo.boundingbox.Diam();
2492 if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
2493 mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
2495 // Local size on faces
2496 occgeo.face_maxh = mparams.maxh;
2498 // Let netgen create _ngMesh and calculate element size on not meshed shapes
2502 int startWith = netgen::MESHCONST_ANALYSE;
2503 int endWith = netgen::MESHCONST_ANALYSE;
2508 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2510 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2512 if(netgen::multithread.terminate)
2515 comment << text(err);
2517 catch (Standard_Failure& ex)
2519 comment << text(ex);
2521 err = 0; //- MESHCONST_ANALYSE isn't so important step
2524 ngLib.setMesh(( Ng_Mesh*) _ngMesh );
2526 _ngMesh->ClearFaceDescriptors(); // we make descriptors our-self
2530 // Pass 1D simple parameters to NETGEN
2531 // --------------------------------
2532 int nbSeg = _simpleHyp->GetNumberOfSegments();
2533 double segSize = _simpleHyp->GetLocalLength();
2534 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
2536 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
2538 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
2539 setLocalSize( e, segSize, *_ngMesh );
2542 else // if ( ! _simpleHyp )
2544 // Local size on vertices and edges
2545 // --------------------------------
2546 for(std::map<int,double>::const_iterator it=EdgeId2LocalSize.begin(); it!=EdgeId2LocalSize.end(); it++)
2548 int key = (*it).first;
2549 double hi = (*it).second;
2550 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2551 const TopoDS_Edge& e = TopoDS::Edge(shape);
2552 setLocalSize( e, hi, *_ngMesh );
2554 for(std::map<int,double>::const_iterator it=VertexId2LocalSize.begin(); it!=VertexId2LocalSize.end(); it++)
2556 int key = (*it).first;
2557 double hi = (*it).second;
2558 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2559 const TopoDS_Vertex& v = TopoDS::Vertex(shape);
2560 gp_Pnt p = BRep_Tool::Pnt(v);
2561 NETGENPlugin_Mesher::RestrictLocalSize( *_ngMesh, p.XYZ(), hi );
2563 for(map<int,double>::const_iterator it=FaceId2LocalSize.begin(); it!=FaceId2LocalSize.end(); it++)
2565 int key = (*it).first;
2566 double val = (*it).second;
2567 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2568 int faceNgID = occgeo.fmap.FindIndex(shape);
2569 if ( faceNgID >= 1 )
2571 occgeo.SetFaceMaxH(faceNgID, val);
2572 for ( TopExp_Explorer edgeExp( shape, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
2573 setLocalSize( TopoDS::Edge( edgeExp.Current() ), val, *_ngMesh );
2577 std::vector<SMESHUtils::ControlPnt> pnt;
2578 SMESHUtils::createPointsSampleFromFace( TopoDS::Face( shape ), val, pnt );
2580 NETGENPlugin_Mesher::RestrictLocalSize( *_ngMesh, pnt[0].XYZ(), val );
2581 for ( size_t i = 1; i < pnt.size(); ++i )
2582 _ngMesh->RestrictLocalH( netgen::Point3d( pnt[i].X(), pnt[i].Y(), pnt[i].Z() ), val );
2585 for(map<int,double>::const_iterator it=SolidId2LocalSize.begin(); it!=SolidId2LocalSize.end(); it++)
2587 int key = (*it).first;
2588 double val = (*it).second;
2589 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2591 std::vector<SMESHUtils::ControlPnt> pnt;
2592 SMESHUtils::createPointsSampleFromSolid( TopoDS::Solid( shape ), val, pnt );
2594 NETGENPlugin_Mesher::RestrictLocalSize( *_ngMesh, pnt[0].XYZ(), val );
2595 for ( size_t i = 1; i < pnt.size(); ++i )
2596 _ngMesh->RestrictLocalH( netgen::Point3d( pnt[i].X(), pnt[i].Y(), pnt[i].Z() ), val );
2600 // Precompute internal edges (issue 0020676) in order to
2601 // add mesh on them correctly (twice) to netgen mesh
2602 if ( !err && internals.hasInternalEdges() )
2604 // load internal shapes into OCCGeometry
2605 netgen::OCCGeometry intOccgeo;
2606 internals.getInternalEdges( intOccgeo.fmap, intOccgeo.emap, intOccgeo.vmap, meshedSM );
2607 intOccgeo.boundingbox = occgeo.boundingbox;
2608 intOccgeo.shape = occgeo.shape;
2609 intOccgeo.face_maxh.SetSize(intOccgeo.fmap.Extent());
2610 intOccgeo.face_maxh = netgen::mparam.maxh;
2611 netgen::Mesh *tmpNgMesh = NULL;
2615 // compute local H on internal shapes in the main mesh
2616 //OCCSetLocalMeshSize(intOccgeo, *_ngMesh); it deletes _ngMesh->localH
2618 // let netgen create a temporary mesh
2620 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2622 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2624 if(netgen::multithread.terminate)
2627 // copy LocalH from the main to temporary mesh
2628 initState.transferLocalH( _ngMesh, tmpNgMesh );
2630 // compute mesh on internal edges
2631 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2633 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2635 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2637 comment << text(err);
2639 catch (Standard_Failure& ex)
2641 comment << text(ex);
2644 initState.restoreLocalH( tmpNgMesh );
2646 // fill SMESH by netgen mesh
2647 vector< const SMDS_MeshNode* > tmpNodeVec;
2648 FillSMesh( intOccgeo, *tmpNgMesh, initState, *_mesh, tmpNodeVec, comment );
2649 err = ( err || !comment.empty() );
2651 nglib::Ng_DeleteMesh((nglib::Ng_Mesh*)tmpNgMesh);
2654 // Fill _ngMesh with nodes and segments of computed submeshes
2657 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_0D ]) &&
2658 FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_1D ], &quadHelper));
2660 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2665 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2670 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2672 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2674 if(netgen::multithread.terminate)
2677 comment << text(err);
2679 catch (Standard_Failure& ex)
2681 comment << text(ex);
2686 _ticTime = ( doneTime += edgeMeshingTime ) / _totalTime / _progressTic;
2688 mparams.uselocalh = true; // restore as it is used at surface optimization
2690 // ---------------------
2691 // compute surface mesh
2692 // ---------------------
2695 // Pass 2D simple parameters to NETGEN
2697 if ( double area = _simpleHyp->GetMaxElementArea() ) {
2699 mparams.maxh = sqrt(2. * area/sqrt(3.0));
2700 mparams.grading = 0.4; // moderate size growth
2703 // length from edges
2704 if ( _ngMesh->GetNSeg() ) {
2705 double edgeLength = 0;
2706 TopTools_MapOfShape visitedEdges;
2707 for ( TopExp_Explorer exp( _shape, TopAbs_EDGE ); exp.More(); exp.Next() )
2708 if( visitedEdges.Add(exp.Current()) )
2709 edgeLength += SMESH_Algo::EdgeLength( TopoDS::Edge( exp.Current() ));
2710 // we have to multiply length by 2 since for each TopoDS_Edge there
2711 // are double set of NETGEN edges, in other words, we have to
2712 // divide _ngMesh->GetNSeg() by 2.
2713 mparams.maxh = 2*edgeLength / _ngMesh->GetNSeg();
2716 mparams.maxh = 1000;
2718 mparams.grading = 0.2; // slow size growth
2720 mparams.quad = _simpleHyp->GetAllowQuadrangles();
2721 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2722 _ngMesh->SetGlobalH (mparams.maxh);
2723 netgen::Box<3> bb = occgeo.GetBoundingBox();
2724 bb.Increase (bb.Diam()/20);
2725 _ngMesh->SetLocalH (bb.PMin(), bb.PMax(), mparams.grading);
2728 // Care of vertices internal in faces (issue 0020676)
2729 if ( internals.hasInternalVertexInFace() )
2731 // store computed segments in SMESH in order not to create SMESH
2732 // edges for ng segments added by AddIntVerticesInFaces()
2733 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2734 // add segments to faces with internal vertices
2735 AddIntVerticesInFaces( occgeo, *_ngMesh, nodeVec, internals );
2736 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2739 // Build viscous layers
2740 if ( _isViscousLayers2D ||
2741 StdMeshers_ViscousLayers2D::HasProxyMesh( TopoDS::Face( occgeo.fmap(1) ), *_mesh ))
2743 if ( !internals.hasInternalVertexInFace() ) {
2744 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2745 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2747 SMESH_ProxyMesh::Ptr viscousMesh;
2748 SMESH_MesherHelper helper( *_mesh );
2749 for ( int faceID = 1; faceID <= occgeo.fmap.Extent(); ++faceID )
2751 const TopoDS_Face& F = TopoDS::Face( occgeo.fmap( faceID ));
2752 viscousMesh = StdMeshers_ViscousLayers2D::Compute( *_mesh, F );
2755 if ( viscousMesh->NbProxySubMeshes() == 0 )
2757 // exclude from computation ng segments built on EDGEs of F
2758 for (int i = 1; i <= _ngMesh->GetNSeg(); i++)
2760 netgen::Segment & seg = _ngMesh->LineSegment(i);
2761 if (seg.si == faceID)
2764 // add new segments to _ngMesh instead of excluded ones
2765 helper.SetSubShape( F );
2767 StdMeshers_FaceSide::GetFaceWires( F, *_mesh, /*skipMediumNodes=*/true,
2768 error, viscousMesh );
2769 error = AddSegmentsToMesh( *_ngMesh, occgeo, wires, helper, nodeVec );
2771 if ( !error ) error = SMESH_ComputeError::New();
2773 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2776 // Let netgen compute 2D mesh
2777 startWith = netgen::MESHCONST_MESHSURFACE;
2778 endWith = _optimize ? netgen::MESHCONST_OPTSURFACE : netgen::MESHCONST_MESHSURFACE;
2783 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2785 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2787 if(netgen::multithread.terminate)
2790 comment << text (err);
2792 catch (Standard_Failure& ex)
2794 comment << text(ex);
2795 //err = 1; -- try to make volumes anyway
2797 catch (netgen::NgException exc)
2799 comment << text(exc);
2800 //err = 1; -- try to make volumes anyway
2805 doneTime += faceMeshingTime + ( _optimize ? faceOptimizTime : 0 );
2806 _ticTime = doneTime / _totalTime / _progressTic;
2808 // ---------------------
2809 // generate volume mesh
2810 // ---------------------
2811 // Fill _ngMesh with nodes and faces of computed 2D submeshes
2812 if ( !err && _isVolume && ( !meshedSM[ MeshDim_2D ].empty() || mparams.quad ))
2814 // load SMESH with computed segments and faces
2815 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2817 // compute pyramids on quadrangles
2818 SMESH_ProxyMesh::Ptr proxyMesh;
2819 if ( _mesh->NbQuadrangles() > 0 )
2820 for ( int iS = 1; iS <= occgeo.somap.Extent(); ++iS )
2822 StdMeshers_QuadToTriaAdaptor* Adaptor = new StdMeshers_QuadToTriaAdaptor;
2823 proxyMesh.reset( Adaptor );
2825 int nbPyrams = _mesh->NbPyramids();
2826 Adaptor->Compute( *_mesh, occgeo.somap(iS) );
2827 if ( nbPyrams != _mesh->NbPyramids() )
2829 list< SMESH_subMesh* > quadFaceSM;
2830 for (TopExp_Explorer face(occgeo.somap(iS), TopAbs_FACE); face.More(); face.Next())
2831 if ( Adaptor->GetProxySubMesh( face.Current() ))
2833 quadFaceSM.push_back( _mesh->GetSubMesh( face.Current() ));
2834 meshedSM[ MeshDim_2D ].remove( quadFaceSM.back() );
2836 FillNgMesh(occgeo, *_ngMesh, nodeVec, quadFaceSM, &quadHelper, proxyMesh);
2839 // fill _ngMesh with faces of sub-meshes
2840 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_2D ], &quadHelper));
2841 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2842 //toPython( _ngMesh, "/tmp/ngPython.py");
2844 if (!err && _isVolume)
2846 // Pass 3D simple parameters to NETGEN
2847 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
2848 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
2850 if ( double vol = simple3d->GetMaxElementVolume() ) {
2852 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
2853 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2856 // length from faces
2857 mparams.maxh = _ngMesh->AverageH();
2859 _ngMesh->SetGlobalH (mparams.maxh);
2860 mparams.grading = 0.4;
2862 _ngMesh->CalcLocalH(mparams.grading);
2864 _ngMesh->CalcLocalH();
2867 // Care of vertices internal in solids and internal faces (issue 0020676)
2868 if ( internals.hasInternalVertexInSolid() || internals.hasInternalFaces() )
2870 // store computed faces in SMESH in order not to create SMESH
2871 // faces for ng faces added here
2872 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2873 // add ng faces to solids with internal vertices
2874 AddIntVerticesInSolids( occgeo, *_ngMesh, nodeVec, internals );
2875 // duplicate mesh faces on internal faces
2876 FixIntFaces( occgeo, *_ngMesh, internals );
2877 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2879 // Let netgen compute 3D mesh
2880 startWith = endWith = netgen::MESHCONST_MESHVOLUME;
2885 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2887 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2889 if(netgen::multithread.terminate)
2892 if ( comment.empty() ) // do not overwrite a previos error
2893 comment << text(err);
2895 catch (Standard_Failure& ex)
2897 if ( comment.empty() ) // do not overwrite a previos error
2898 comment << text(ex);
2901 catch (netgen::NgException exc)
2903 if ( comment.empty() ) // do not overwrite a previos error
2904 comment << text(exc);
2907 _ticTime = ( doneTime += voluMeshingTime ) / _totalTime / _progressTic;
2909 // Let netgen optimize 3D mesh
2910 if ( !err && _optimize )
2912 startWith = endWith = netgen::MESHCONST_OPTVOLUME;
2917 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2919 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2921 if(netgen::multithread.terminate)
2924 if ( comment.empty() ) // do not overwrite a previos error
2925 comment << text(err);
2927 catch (Standard_Failure& ex)
2929 if ( comment.empty() ) // do not overwrite a previos error
2930 comment << text(ex);
2932 catch (netgen::NgException exc)
2934 if ( comment.empty() ) // do not overwrite a previos error
2935 comment << text(exc);
2939 if (!err && mparams.secondorder > 0)
2944 if ( !meshedSM[ MeshDim_1D ].empty() )
2946 // remove segments not attached to geometry (IPAL0052479)
2947 for (int i = 1; i <= _ngMesh->GetNSeg(); ++i)
2949 const netgen::Segment & seg = _ngMesh->LineSegment (i);
2950 if ( seg.epgeominfo[ 0 ].edgenr == 0 )
2951 _ngMesh->DeleteSegment( i );
2953 _ngMesh->Compress();
2955 // convert to quadratic
2956 netgen::OCCRefinementSurfaces ref (occgeo);
2957 ref.MakeSecondOrder (*_ngMesh);
2959 // care of elements already loaded to SMESH
2960 // if ( initState._nbSegments > 0 )
2961 // makeQuadratic( occgeo.emap, _mesh );
2962 // if ( initState._nbFaces > 0 )
2963 // makeQuadratic( occgeo.fmap, _mesh );
2965 catch (Standard_Failure& ex)
2967 if ( comment.empty() ) // do not overwrite a previos error
2968 comment << "Exception in netgen at passing to 2nd order ";
2970 catch (netgen::NgException exc)
2972 if ( comment.empty() ) // do not overwrite a previos error
2973 comment << exc.What();
2978 _ticTime = 0.98 / _progressTic;
2980 //int nbNod = _ngMesh->GetNP();
2981 //int nbSeg = _ngMesh->GetNSeg();
2982 int nbFac = _ngMesh->GetNSE();
2983 int nbVol = _ngMesh->GetNE();
2984 bool isOK = ( !err && (_isVolume ? (nbVol > 0) : (nbFac > 0)) );
2986 // Feed back the SMESHDS with the generated Nodes and Elements
2987 if ( true /*isOK*/ ) // get whatever built
2989 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2991 if ( quadHelper.GetIsQuadratic() ) // remove free nodes
2992 for ( size_t i = 0; i < nodeVec.size(); ++i )
2993 if ( nodeVec[i] && nodeVec[i]->NbInverseElements() == 0 )
2994 _mesh->GetMeshDS()->RemoveFreeNode( nodeVec[i], 0, /*fromGroups=*/false );
2996 SMESH_ComputeErrorPtr readErr = ReadErrors(nodeVec);
2997 if ( readErr && !readErr->myBadElements.empty() )
3000 if ( !comment.empty() && !readErr->myComment.empty() ) comment += "\n";
3001 comment += readErr->myComment;
3003 if ( error->IsOK() && ( !isOK || comment.size() > 0 ))
3004 error->myName = COMPERR_ALGO_FAILED;
3005 if ( !comment.empty() )
3006 error->myComment = comment;
3008 // SetIsAlwaysComputed( true ) to empty sub-meshes, which
3009 // appear if the geometry contains coincident sub-shape due
3010 // to bool merge_solids = 1; in netgen/libsrc/occ/occgenmesh.cpp
3011 const int nbMaps = 2;
3012 const TopTools_IndexedMapOfShape* geoMaps[nbMaps] =
3013 { & occgeo.vmap, & occgeo.emap/*, & occgeo.fmap*/ };
3014 for ( int iMap = 0; iMap < nbMaps; ++iMap )
3015 for (int i = 1; i <= geoMaps[iMap]->Extent(); i++)
3016 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( geoMaps[iMap]->FindKey(i)))
3017 if ( !sm->IsMeshComputed() )
3018 sm->SetIsAlwaysComputed( true );
3020 // set bad compute error to subshapes of all failed sub-shapes
3021 if ( !error->IsOK() )
3023 bool pb2D = false, pb3D = false;
3024 for (int i = 1; i <= occgeo.fmap.Extent(); i++) {
3025 int status = occgeo.facemeshstatus[i-1];
3026 if (status == netgen::FACE_MESHED_OK ) continue;
3027 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.fmap( i ))) {
3028 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
3029 if ( !smError || smError->IsOK() ) {
3030 if ( status == netgen::FACE_FAILED )
3031 smError.reset( new SMESH_ComputeError( *error ));
3033 smError.reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED, "Ignored" ));
3034 if ( SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
3035 smError->myName = COMPERR_WARNING;
3037 pb2D = pb2D || smError->IsKO();
3040 if ( !pb2D ) // all faces are OK
3041 for (int i = 1; i <= occgeo.somap.Extent(); i++)
3042 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.somap( i )))
3044 bool smComputed = nbVol && !sm->IsEmpty();
3045 if ( smComputed && internals.hasInternalVertexInSolid( sm->GetId() ))
3047 int nbIntV = internals.getSolidsWithVertices().find( sm->GetId() )->second.size();
3048 SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
3049 smComputed = ( smDS->NbElements() > 0 || smDS->NbNodes() > nbIntV );
3051 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
3052 if ( !smComputed && ( !smError || smError->IsOK() ))
3054 smError.reset( new SMESH_ComputeError( *error ));
3055 if ( nbVol && SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
3057 smError->myName = COMPERR_WARNING;
3059 else if ( !smError->myBadElements.empty() ) // bad surface mesh
3061 if ( !hasBadElemOnSolid( smError->myBadElements, sm ))
3065 pb3D = pb3D || ( smError && smError->IsKO() );
3067 if ( !pb2D && !pb3D )
3068 err = 0; // no fatal errors, only warnings
3071 ngLib._isComputeOk = !err;
3076 //=============================================================================
3080 //=============================================================================
3081 bool NETGENPlugin_Mesher::Evaluate(MapShapeNbElems& aResMap)
3083 netgen::MeshingParameters& mparams = netgen::mparam;
3086 // -------------------------
3087 // Prepare OCC geometry
3088 // -------------------------
3089 netgen::OCCGeometry occgeo;
3090 list< SMESH_subMesh* > meshedSM[4]; // for 0-3 dimensions
3091 NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
3092 PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
3094 bool tooManyElems = false;
3095 const int hugeNb = std::numeric_limits<int>::max() / 100;
3100 // pass 1D simple parameters to NETGEN
3103 // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
3104 mparams.uselocalh = false;
3105 mparams.grading = 0.8; // not limitited size growth
3107 if ( _simpleHyp->GetNumberOfSegments() )
3109 mparams.maxh = occgeo.boundingbox.Diam();
3112 mparams.maxh = _simpleHyp->GetLocalLength();
3115 if ( mparams.maxh == 0.0 )
3116 mparams.maxh = occgeo.boundingbox.Diam();
3117 if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
3118 mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
3120 // let netgen create _ngMesh and calculate element size on not meshed shapes
3121 NETGENPlugin_NetgenLibWrapper ngLib;
3122 netgen::Mesh *ngMesh = NULL;
3126 int startWith = netgen::MESHCONST_ANALYSE;
3127 int endWith = netgen::MESHCONST_MESHEDGES;
3129 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, mparams, startWith, endWith);
3131 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, startWith, endWith, optstr);
3134 if(netgen::multithread.terminate)
3137 ngLib.setMesh(( Ng_Mesh*) ngMesh );
3139 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( _shape ))
3140 sm->GetComputeError().reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED ));
3145 // Pass 1D simple parameters to NETGEN
3146 // --------------------------------
3147 int nbSeg = _simpleHyp->GetNumberOfSegments();
3148 double segSize = _simpleHyp->GetLocalLength();
3149 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
3151 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
3153 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
3154 setLocalSize( e, segSize, *ngMesh );
3157 else // if ( ! _simpleHyp )
3159 // Local size on vertices and edges
3160 // --------------------------------
3161 for(std::map<int,double>::const_iterator it=EdgeId2LocalSize.begin(); it!=EdgeId2LocalSize.end(); it++)
3163 int key = (*it).first;
3164 double hi = (*it).second;
3165 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3166 const TopoDS_Edge& e = TopoDS::Edge(shape);
3167 setLocalSize( e, hi, *ngMesh );
3169 for(std::map<int,double>::const_iterator it=VertexId2LocalSize.begin(); it!=VertexId2LocalSize.end(); it++)
3171 int key = (*it).first;
3172 double hi = (*it).second;
3173 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3174 const TopoDS_Vertex& v = TopoDS::Vertex(shape);
3175 gp_Pnt p = BRep_Tool::Pnt(v);
3176 NETGENPlugin_Mesher::RestrictLocalSize( *ngMesh, p.XYZ(), hi );
3178 for(map<int,double>::const_iterator it=FaceId2LocalSize.begin();
3179 it!=FaceId2LocalSize.end(); it++)
3181 int key = (*it).first;
3182 double val = (*it).second;
3183 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3184 int faceNgID = occgeo.fmap.FindIndex(shape);
3185 occgeo.SetFaceMaxH(faceNgID, val);
3186 for ( TopExp_Explorer edgeExp( shape, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
3187 setLocalSize( TopoDS::Edge( edgeExp.Current() ), val, *ngMesh );
3190 // calculate total nb of segments and length of edges
3191 double fullLen = 0.0;
3193 int entity = mparams.secondorder > 0 ? SMDSEntity_Quad_Edge : SMDSEntity_Edge;
3194 TopTools_DataMapOfShapeInteger Edge2NbSeg;
3195 for (TopExp_Explorer exp(_shape, TopAbs_EDGE); exp.More(); exp.Next())
3197 TopoDS_Edge E = TopoDS::Edge( exp.Current() );
3198 if( !Edge2NbSeg.Bind(E,0) )
3201 double aLen = SMESH_Algo::EdgeLength(E);
3204 vector<int>& aVec = aResMap[_mesh->GetSubMesh(E)];
3206 aVec.resize( SMDSEntity_Last, 0);
3208 fullNbSeg += aVec[ entity ];
3211 // store nb of segments computed by Netgen
3212 NCollection_Map<Link> linkMap;
3213 for (int i = 1; i <= ngMesh->GetNSeg(); ++i )
3215 const netgen::Segment& seg = ngMesh->LineSegment(i);
3216 Link link(seg[0], seg[1]);
3217 if ( !linkMap.Add( link )) continue;
3218 int aGeomEdgeInd = seg.epgeominfo[0].edgenr;
3219 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
3221 vector<int>& aVec = aResMap[_mesh->GetSubMesh(occgeo.emap(aGeomEdgeInd))];
3225 // store nb of nodes on edges computed by Netgen
3226 TopTools_DataMapIteratorOfDataMapOfShapeInteger Edge2NbSegIt(Edge2NbSeg);
3227 for (; Edge2NbSegIt.More(); Edge2NbSegIt.Next())
3229 vector<int>& aVec = aResMap[_mesh->GetSubMesh(Edge2NbSegIt.Key())];
3230 if ( aVec[ entity ] > 1 && aVec[ SMDSEntity_Node ] == 0 )
3231 aVec[SMDSEntity_Node] = mparams.secondorder > 0 ? 2*aVec[ entity ]-1 : aVec[ entity ]-1;
3233 fullNbSeg += aVec[ entity ];
3234 Edge2NbSeg( Edge2NbSegIt.Key() ) = aVec[ entity ];
3236 if ( fullNbSeg == 0 )
3243 if ( double area = _simpleHyp->GetMaxElementArea() ) {
3245 mparams.maxh = sqrt(2. * area/sqrt(3.0));
3246 mparams.grading = 0.4; // moderate size growth
3249 // length from edges
3250 mparams.maxh = fullLen/fullNbSeg;
3251 mparams.grading = 0.2; // slow size growth
3254 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3255 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3257 for (TopExp_Explorer exp(_shape, TopAbs_FACE); exp.More(); exp.Next())
3259 TopoDS_Face F = TopoDS::Face( exp.Current() );
3260 SMESH_subMesh *sm = _mesh->GetSubMesh(F);
3262 BRepGProp::SurfaceProperties(F,G);
3263 double anArea = G.Mass();
3264 tooManyElems = tooManyElems || ( anArea/hugeNb > mparams.maxh*mparams.maxh );
3266 if ( !tooManyElems )
3268 TopTools_MapOfShape egdes;
3269 for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next())
3270 if ( egdes.Add( exp1.Current() ))
3271 nb1d += Edge2NbSeg.Find(exp1.Current());
3273 int nbFaces = tooManyElems ? hugeNb : int( 4*anArea / (mparams.maxh*mparams.maxh*sqrt(3.)));
3274 int nbNodes = tooManyElems ? hugeNb : (( nbFaces*3 - (nb1d-1)*2 ) / 6 + 1 );
3276 vector<int> aVec(SMDSEntity_Last, 0);
3277 if( mparams.secondorder > 0 ) {
3278 int nb1d_in = (nbFaces*3 - nb1d) / 2;
3279 aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
3280 aVec[SMDSEntity_Quad_Triangle] = nbFaces;
3283 aVec[SMDSEntity_Node] = Max ( nbNodes, 0 );
3284 aVec[SMDSEntity_Triangle] = nbFaces;
3286 aResMap[sm].swap(aVec);
3293 // pass 3D simple parameters to NETGEN
3294 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
3295 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
3297 if ( double vol = simple3d->GetMaxElementVolume() ) {
3299 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
3300 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3303 // using previous length from faces
3305 mparams.grading = 0.4;
3306 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3309 BRepGProp::VolumeProperties(_shape,G);
3310 double aVolume = G.Mass();
3311 double tetrVol = 0.1179*mparams.maxh*mparams.maxh*mparams.maxh;
3312 tooManyElems = tooManyElems || ( aVolume/hugeNb > tetrVol );
3313 int nbVols = tooManyElems ? hugeNb : int(aVolume/tetrVol);
3314 int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
3315 vector<int> aVec(SMDSEntity_Last, 0 );
3316 if ( tooManyElems ) // avoid FPE
3318 aVec[SMDSEntity_Node] = hugeNb;
3319 aVec[ mparams.secondorder > 0 ? SMDSEntity_Quad_Tetra : SMDSEntity_Tetra] = hugeNb;
3323 if( mparams.secondorder > 0 ) {
3324 aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
3325 aVec[SMDSEntity_Quad_Tetra] = nbVols;
3328 aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
3329 aVec[SMDSEntity_Tetra] = nbVols;
3332 SMESH_subMesh *sm = _mesh->GetSubMesh(_shape);
3333 aResMap[sm].swap(aVec);
3339 double NETGENPlugin_Mesher::GetProgress(const SMESH_Algo* holder,
3340 const int * algoProgressTic,
3341 const double * algoProgress) const
3343 ((int&) _progressTic ) = *algoProgressTic + 1;
3345 if ( !_occgeom ) return 0;
3347 double progress = -1;
3350 if ( _ticTime < 0 && netgen::multithread.task[0] == 'O'/*Optimizing surface*/ )
3352 ((double&) _ticTime ) = edgeFaceMeshingTime / _totalTime / _progressTic;
3354 else if ( !_optimize /*&& _occgeom->fmap.Extent() > 1*/ )
3356 int doneShapeIndex = -1;
3357 while ( doneShapeIndex+1 < _occgeom->facemeshstatus.Size() &&
3358 _occgeom->facemeshstatus[ doneShapeIndex+1 ])
3360 if ( doneShapeIndex+1 != _curShapeIndex )
3362 ((int&) _curShapeIndex) = doneShapeIndex+1;
3363 double doneShapeRate = _curShapeIndex / double( _occgeom->fmap.Extent() );
3364 double doneTime = edgeMeshingTime + doneShapeRate * faceMeshingTime;
3365 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3366 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3367 // << " " << doneTime / _totalTime / _progressTic << endl;
3371 else if ( !_optimize && _occgeom->somap.Extent() > 1 )
3373 int curShapeIndex = _curShapeIndex;
3374 if ( _ngMesh->GetNE() > 0 )
3376 netgen::Element el = (*_ngMesh)[netgen::ElementIndex( _ngMesh->GetNE()-1 )];
3377 curShapeIndex = el.GetIndex();
3379 if ( curShapeIndex != _curShapeIndex )
3381 ((int&) _curShapeIndex) = curShapeIndex;
3382 double doneShapeRate = _curShapeIndex / double( _occgeom->somap.Extent() );
3383 double doneTime = edgeFaceMeshingTime + doneShapeRate * voluMeshingTime;
3384 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3385 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3386 // << " " << doneTime / _totalTime / _progressTic << endl;
3390 progress = Max( *algoProgressTic * _ticTime, *algoProgress );
3393 ((int&) *algoProgressTic )++;
3394 ((double&) *algoProgress) = progress;
3396 //cout << progress << " " << *algoProgressTic << " " << netgen::multithread.task << " "<< _ticTime << endl;
3398 return Min( progress, 0.99 );
3401 //================================================================================
3403 * \brief Read mesh entities preventing successful computation from "test.out" file
3405 //================================================================================
3407 SMESH_ComputeErrorPtr
3408 NETGENPlugin_Mesher::ReadErrors(const vector<const SMDS_MeshNode* >& nodeVec)
3410 SMESH_ComputeErrorPtr err = SMESH_ComputeError::New
3411 (COMPERR_BAD_INPUT_MESH, "Some edges multiple times in surface mesh");
3412 SMESH_File file("test.out");
3414 vector<int> three1(3), three2(3);
3415 const char* badEdgeStr = " multiple times in surface mesh";
3416 const int badEdgeStrLen = strlen( badEdgeStr );
3417 const int nbNodes = nodeVec.size();
3419 while( !file.eof() )
3421 if ( strncmp( file, "Edge ", 5 ) == 0 &&
3422 file.getInts( two ) &&
3423 strncmp( file, badEdgeStr, badEdgeStrLen ) == 0 &&
3424 two[0] < nbNodes && two[1] < nbNodes )
3426 err->myBadElements.push_back( new SMDS_LinearEdge( nodeVec[ two[0]], nodeVec[ two[1]] ));
3427 file += badEdgeStrLen;
3429 else if ( strncmp( file, "Intersecting: ", 14 ) == 0 )
3432 // openelement 18 with open element 126
3436 const char* pos = file;
3437 bool ok = ( strncmp( file, "openelement ", 12 ) == 0 );
3438 ok = ok && file.getInts( two );
3439 ok = ok && file.getInts( three1 );
3440 ok = ok && file.getInts( three2 );
3441 for ( int i = 0; ok && i < 3; ++i )
3442 ok = ( three1[i] < nbNodes && nodeVec[ three1[i]]);
3443 for ( int i = 0; ok && i < 3; ++i )
3444 ok = ( three2[i] < nbNodes && nodeVec[ three2[i]]);
3447 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three1[0]],
3448 nodeVec[ three1[1]],
3449 nodeVec[ three1[2]]));
3450 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three2[0]],
3451 nodeVec[ three2[1]],
3452 nodeVec[ three2[2]]));
3453 err->myComment = "Intersecting triangles";
3467 size_t nbBadElems = err->myBadElements.size();
3468 if ( nbBadElems ) nbBadElems++; // avoid warning: variable set but not used
3474 //================================================================================
3476 * \brief Write a python script creating an equivalent SALOME mesh.
3477 * This is useful to see what mesh is passed as input for the next step of mesh
3478 * generation (of mesh of higher dimension)
3480 //================================================================================
3482 void NETGENPlugin_Mesher::toPython( const netgen::Mesh* ngMesh,
3483 const std::string& pyFile)
3485 ofstream outfile(pyFile.c_str(), ios::out);
3486 if ( !outfile ) return;
3488 outfile << "import SMESH" << endl
3489 << "from salome.smesh import smeshBuilder" << endl
3490 << "smesh = smeshBuilder.New(salome.myStudy)" << endl
3491 << "mesh = smesh.Mesh()" << endl << endl;
3493 using namespace netgen;
3495 for (pi = PointIndex::BASE;
3496 pi < ngMesh->GetNP()+PointIndex::BASE; pi++)
3498 outfile << "mesh.AddNode( ";
3499 outfile << (*ngMesh)[pi](0) << ", ";
3500 outfile << (*ngMesh)[pi](1) << ", ";
3501 outfile << (*ngMesh)[pi](2) << ") ## "<< pi << endl;
3504 int nbDom = ngMesh->GetNDomains();
3505 for ( int i = 0; i < nbDom; ++i )
3506 outfile<< "grp" << i+1 << " = mesh.CreateEmptyGroup( SMESH.FACE, 'domain"<< i+1 << "')"<< endl;
3508 SurfaceElementIndex sei;
3509 for (sei = 0; sei < ngMesh->GetNSE(); sei++)
3511 outfile << "mesh.AddFace([ ";
3512 Element2d sel = (*ngMesh)[sei];
3513 for (int j = 0; j < sel.GetNP(); j++)
3514 outfile << sel[j] << ( j+1 < sel.GetNP() ? ", " : " ])");
3515 if ( sel.IsDeleted() ) outfile << " ## IsDeleted ";
3518 if ((*ngMesh)[sei].GetIndex())
3520 if ( int dom1 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainIn())
3521 outfile << "grp"<< dom1 <<".Add([ " << (int)sei+1 << " ])" << endl;
3522 if ( int dom2 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainOut())
3523 outfile << "grp"<< dom2 <<".Add([ " << (int)sei+1 << " ])" << endl;
3527 for (ElementIndex ei = 0; ei < ngMesh->GetNE(); ei++)
3529 Element el = (*ngMesh)[ei];
3530 outfile << "mesh.AddVolume([ ";
3531 for (int j = 0; j < el.GetNP(); j++)
3532 outfile << el[j] << ( j+1 < el.GetNP() ? ", " : " ])");
3536 for (int i = 1; i <= ngMesh->GetNSeg(); i++)
3538 const Segment & seg = ngMesh->LineSegment (i);
3539 outfile << "mesh.AddEdge([ "
3541 << seg[1] << " ])" << endl;
3543 cout << "Write " << pyFile << endl;
3546 //================================================================================
3548 * \brief Constructor of NETGENPlugin_ngMeshInfo
3550 //================================================================================
3552 NETGENPlugin_ngMeshInfo::NETGENPlugin_ngMeshInfo( netgen::Mesh* ngMesh):
3557 _nbNodes = ngMesh->GetNP();
3558 _nbSegments = ngMesh->GetNSeg();
3559 _nbFaces = ngMesh->GetNSE();
3560 _nbVolumes = ngMesh->GetNE();
3564 _nbNodes = _nbSegments = _nbFaces = _nbVolumes = 0;
3568 //================================================================================
3570 * \brief Copy LocalH member from one netgen mesh to another
3572 //================================================================================
3574 void NETGENPlugin_ngMeshInfo::transferLocalH( netgen::Mesh* fromMesh,
3575 netgen::Mesh* toMesh )
3577 if ( !fromMesh->LocalHFunctionGenerated() ) return;
3578 if ( !toMesh->LocalHFunctionGenerated() )
3580 toMesh->CalcLocalH(netgen::mparam.grading);
3582 toMesh->CalcLocalH();
3585 const size_t size = sizeof( netgen::LocalH );
3586 _copyOfLocalH = new char[ size ];
3587 memcpy( (void*)_copyOfLocalH, (void*)&toMesh->LocalHFunction(), size );
3588 memcpy( (void*)&toMesh->LocalHFunction(), (void*)&fromMesh->LocalHFunction(), size );
3591 //================================================================================
3593 * \brief Restore LocalH member of a netgen mesh
3595 //================================================================================
3597 void NETGENPlugin_ngMeshInfo::restoreLocalH( netgen::Mesh* toMesh )
3599 if ( _copyOfLocalH )
3601 const size_t size = sizeof( netgen::LocalH );
3602 memcpy( (void*)&toMesh->LocalHFunction(), (void*)_copyOfLocalH, size );
3603 delete [] _copyOfLocalH;
3608 //================================================================================
3610 * \brief Find "internal" sub-shapes
3612 //================================================================================
3614 NETGENPlugin_Internals::NETGENPlugin_Internals( SMESH_Mesh& mesh,
3615 const TopoDS_Shape& shape,
3617 : _mesh( mesh ), _is3D( is3D )
3619 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3621 TopExp_Explorer f,e;
3622 for ( f.Init( shape, TopAbs_FACE ); f.More(); f.Next() )
3624 int faceID = meshDS->ShapeToIndex( f.Current() );
3626 // find not computed internal edges
3628 for ( e.Init( f.Current().Oriented(TopAbs_FORWARD), TopAbs_EDGE ); e.More(); e.Next() )
3629 if ( e.Current().Orientation() == TopAbs_INTERNAL )
3631 SMESH_subMesh* eSM = mesh.GetSubMesh( e.Current() );
3632 if ( eSM->IsEmpty() )
3634 _e2face.insert( make_pair( eSM->GetId(), faceID ));
3635 for ( TopoDS_Iterator v(e.Current()); v.More(); v.Next() )
3636 _e2face.insert( make_pair( meshDS->ShapeToIndex( v.Value() ), faceID ));
3640 // find internal vertices in a face
3641 set<int> intVV; // issue 0020850 where same vertex is twice in a face
3642 for ( TopoDS_Iterator fSub( f.Current() ); fSub.More(); fSub.Next())
3643 if ( fSub.Value().ShapeType() == TopAbs_VERTEX )
3645 int vID = meshDS->ShapeToIndex( fSub.Value() );
3646 if ( intVV.insert( vID ).second )
3647 _f2v[ faceID ].push_back( vID );
3652 // find internal faces and their subshapes where nodes are to be doubled
3653 // to make a crack with non-sewed borders
3655 if ( f.Current().Orientation() == TopAbs_INTERNAL )
3657 _intShapes.insert( meshDS->ShapeToIndex( f.Current() ));
3660 list< TopoDS_Shape > edges;
3661 for ( e.Init( f.Current(), TopAbs_EDGE ); e.More(); e.Next())
3662 if ( SMESH_MesherHelper::NbAncestors( e.Current(), mesh, TopAbs_FACE ) > 1 )
3664 _intShapes.insert( meshDS->ShapeToIndex( e.Current() ));
3665 edges.push_back( e.Current() );
3666 // find border faces
3667 PShapeIteratorPtr fIt =
3668 SMESH_MesherHelper::GetAncestors( edges.back(),mesh,TopAbs_FACE );
3669 while ( const TopoDS_Shape* pFace = fIt->next() )
3670 if ( !pFace->IsSame( f.Current() ))
3671 _borderFaces.insert( meshDS->ShapeToIndex( *pFace ));
3674 // we consider vertex internal if it is shared by more than one internal edge
3675 list< TopoDS_Shape >::iterator edge = edges.begin();
3676 for ( ; edge != edges.end(); ++edge )
3677 for ( TopoDS_Iterator v( *edge ); v.More(); v.Next() )
3679 set<int> internalEdges;
3680 PShapeIteratorPtr eIt =
3681 SMESH_MesherHelper::GetAncestors( v.Value(),mesh,TopAbs_EDGE );
3682 while ( const TopoDS_Shape* pEdge = eIt->next() )
3684 int edgeID = meshDS->ShapeToIndex( *pEdge );
3685 if ( isInternalShape( edgeID ))
3686 internalEdges.insert( edgeID );
3688 if ( internalEdges.size() > 1 )
3689 _intShapes.insert( meshDS->ShapeToIndex( v.Value() ));
3693 } // loop on geom faces
3695 // find vertices internal in solids
3698 for ( TopExp_Explorer so(shape, TopAbs_SOLID); so.More(); so.Next())
3700 int soID = meshDS->ShapeToIndex( so.Current() );
3701 for ( TopoDS_Iterator soSub( so.Current() ); soSub.More(); soSub.Next())
3702 if ( soSub.Value().ShapeType() == TopAbs_VERTEX )
3703 _s2v[ soID ].push_back( meshDS->ShapeToIndex( soSub.Value() ));
3708 //================================================================================
3710 * \brief Find mesh faces on non-internal geom faces sharing internal edge
3711 * some nodes of which are to be doubled to make the second border of the "crack"
3713 //================================================================================
3715 void NETGENPlugin_Internals::findBorderElements( TIDSortedElemSet & borderElems )
3717 if ( _intShapes.empty() ) return;
3719 SMESH_Mesh& mesh = const_cast<SMESH_Mesh&>(_mesh);
3720 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3722 // loop on internal geom edges
3723 set<int>::const_iterator intShapeId = _intShapes.begin();
3724 for ( ; intShapeId != _intShapes.end(); ++intShapeId )
3726 const TopoDS_Shape& s = meshDS->IndexToShape( *intShapeId );
3727 if ( s.ShapeType() != TopAbs_EDGE ) continue;
3729 // get internal and non-internal geom faces sharing the internal edge <s>
3731 set<int>::iterator bordFace = _borderFaces.end();
3732 PShapeIteratorPtr faces = SMESH_MesherHelper::GetAncestors( s, _mesh, TopAbs_FACE );
3733 while ( const TopoDS_Shape* pFace = faces->next() )
3735 int faceID = meshDS->ShapeToIndex( *pFace );
3736 if ( isInternalShape( faceID ))
3739 bordFace = _borderFaces.insert( faceID ).first;
3741 if ( bordFace == _borderFaces.end() || !intFace ) continue;
3743 // get all links of mesh faces on internal geom face sharing nodes on edge <s>
3744 set< SMESH_OrientedLink > links; //!< links of faces on internal geom face
3745 list<const SMDS_MeshElement*> suspectFaces[2]; //!< mesh faces on border geom faces
3746 int nbSuspectFaces = 0;
3747 SMESHDS_SubMesh* intFaceSM = meshDS->MeshElements( intFace );
3748 if ( !intFaceSM || intFaceSM->NbElements() == 0 ) continue;
3749 SMESH_subMeshIteratorPtr smIt = mesh.GetSubMesh( s )->getDependsOnIterator(true,true);
3750 while ( smIt->more() )
3752 SMESHDS_SubMesh* sm = smIt->next()->GetSubMeshDS();
3753 if ( !sm ) continue;
3754 SMDS_NodeIteratorPtr nIt = sm->GetNodes();
3755 while ( nIt->more() )
3757 const SMDS_MeshNode* nOnEdge = nIt->next();
3758 SMDS_ElemIteratorPtr fIt = nOnEdge->GetInverseElementIterator(SMDSAbs_Face);
3759 while ( fIt->more() )
3761 const SMDS_MeshElement* f = fIt->next();
3762 const int nbNodes = f->NbCornerNodes();
3763 if ( intFaceSM->Contains( f ))
3765 for ( int i = 0; i < nbNodes; ++i )
3766 links.insert( SMESH_OrientedLink( f->GetNode(i), f->GetNode((i+1)%nbNodes)));
3771 for ( int i = 0; i < nbNodes; ++i )
3772 nbDblNodes += isInternalShape( f->GetNode(i)->getshapeId() );
3774 suspectFaces[ nbDblNodes < 2 ].push_back( f );
3780 // suspectFaces[0] having link with same orientation as mesh faces on
3781 // the internal geom face are <borderElems>. suspectFaces[1] have
3782 // only one node on edge <s>, we decide on them later (at the 2nd loop)
3783 // by links of <borderElems> found at the 1st and 2nd loops
3784 set< SMESH_OrientedLink > borderLinks;
3785 for ( int isPostponed = 0; isPostponed < 2; ++isPostponed )
3787 list<const SMDS_MeshElement*>::iterator fIt = suspectFaces[isPostponed].begin();
3788 for ( int nbF = 0; fIt != suspectFaces[isPostponed].end(); ++fIt, ++nbF )
3790 const SMDS_MeshElement* f = *fIt;
3791 bool isBorder = false, linkFound = false, borderLinkFound = false;
3792 list< SMESH_OrientedLink > faceLinks;
3793 int nbNodes = f->NbCornerNodes();
3794 for ( int i = 0; i < nbNodes; ++i )
3796 SMESH_OrientedLink link( f->GetNode(i), f->GetNode((i+1)%nbNodes));
3797 faceLinks.push_back( link );
3800 set< SMESH_OrientedLink >::iterator foundLink = links.find( link );
3801 if ( foundLink != links.end() )
3804 isBorder = ( foundLink->_reversed == link._reversed );
3805 if ( !isBorder && !isPostponed ) break;
3806 faceLinks.pop_back();
3808 else if ( isPostponed && !borderLinkFound )
3810 foundLink = borderLinks.find( link );
3811 if ( foundLink != borderLinks.end() )
3813 borderLinkFound = true;
3814 isBorder = ( foundLink->_reversed != link._reversed );
3821 borderElems.insert( f );
3822 borderLinks.insert( faceLinks.begin(), faceLinks.end() );
3824 else if ( !linkFound && !borderLinkFound )
3826 suspectFaces[1].push_back( f );
3827 if ( nbF > 2 * nbSuspectFaces )
3828 break; // dead loop protection
3835 //================================================================================
3837 * \brief put internal shapes in maps and fill in submeshes to precompute
3839 //================================================================================
3841 void NETGENPlugin_Internals::getInternalEdges( TopTools_IndexedMapOfShape& fmap,
3842 TopTools_IndexedMapOfShape& emap,
3843 TopTools_IndexedMapOfShape& vmap,
3844 list< SMESH_subMesh* > smToPrecompute[])
3846 if ( !hasInternalEdges() ) return;
3847 map<int,int>::const_iterator ev_face = _e2face.begin();
3848 for ( ; ev_face != _e2face.end(); ++ev_face )
3850 const TopoDS_Shape& ev = _mesh.GetMeshDS()->IndexToShape( ev_face->first );
3851 const TopoDS_Shape& face = _mesh.GetMeshDS()->IndexToShape( ev_face->second );
3853 ( ev.ShapeType() == TopAbs_EDGE ? emap : vmap ).Add( ev );
3855 //cout<<"INTERNAL EDGE or VERTEX "<<ev_face->first<<" on face "<<ev_face->second<<endl;
3857 smToPrecompute[ MeshDim_1D ].push_back( _mesh.GetSubMeshContaining( ev_face->first ));
3861 //================================================================================
3863 * \brief return shapes and submeshes to be meshed and already meshed boundary submeshes
3865 //================================================================================
3867 void NETGENPlugin_Internals::getInternalFaces( TopTools_IndexedMapOfShape& fmap,
3868 TopTools_IndexedMapOfShape& emap,
3869 list< SMESH_subMesh* >& intFaceSM,
3870 list< SMESH_subMesh* >& boundarySM)
3872 if ( !hasInternalFaces() ) return;
3874 // <fmap> and <emap> are for not yet meshed shapes
3875 // <intFaceSM> is for submeshes of faces
3876 // <boundarySM> is for meshed edges and vertices
3881 set<int> shapeIDs ( _intShapes );
3882 if ( !_borderFaces.empty() )
3883 shapeIDs.insert( _borderFaces.begin(), _borderFaces.end() );
3885 set<int>::const_iterator intS = shapeIDs.begin();
3886 for ( ; intS != shapeIDs.end(); ++intS )
3888 SMESH_subMesh* sm = _mesh.GetSubMeshContaining( *intS );
3890 if ( sm->GetSubShape().ShapeType() != TopAbs_FACE ) continue;
3892 intFaceSM.push_back( sm );
3894 // add submeshes of not computed internal faces
3895 if ( !sm->IsEmpty() ) continue;
3897 SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(true,true);
3898 while ( smIt->more() )
3901 const TopoDS_Shape& s = sm->GetSubShape();
3903 if ( sm->IsEmpty() )
3906 switch ( s.ShapeType() ) {
3907 case TopAbs_FACE: fmap.Add ( s ); break;
3908 case TopAbs_EDGE: emap.Add ( s ); break;
3914 if ( s.ShapeType() != TopAbs_FACE )
3915 boundarySM.push_back( sm );
3921 //================================================================================
3923 * \brief Return true if given shape is to be precomputed in order to be correctly
3924 * added to netgen mesh
3926 //================================================================================
3928 bool NETGENPlugin_Internals::isShapeToPrecompute(const TopoDS_Shape& s)
3930 int shapeID = _mesh.GetMeshDS()->ShapeToIndex( s );
3931 switch ( s.ShapeType() ) {
3932 case TopAbs_FACE : break; //return isInternalShape( shapeID ) || isBorderFace( shapeID );
3933 case TopAbs_EDGE : return isInternalEdge( shapeID );
3934 case TopAbs_VERTEX: break;
3940 //================================================================================
3942 * \brief Return SMESH
3944 //================================================================================
3946 SMESH_Mesh& NETGENPlugin_Internals::getMesh() const
3948 return const_cast<SMESH_Mesh&>( _mesh );
3951 //================================================================================
3953 * \brief Access to a counter of NETGENPlugin_NetgenLibWrapper instances
3955 //================================================================================
3957 int& NETGENPlugin_NetgenLibWrapper::instanceCounter()
3959 static int theCouner = 0;
3963 //================================================================================
3965 * \brief Initialize netgen library
3967 //================================================================================
3969 NETGENPlugin_NetgenLibWrapper::NETGENPlugin_NetgenLibWrapper()
3971 if ( instanceCounter() == 0 )
3974 ++instanceCounter();
3976 _isComputeOk = false;
3980 if ( !getenv( "KEEP_NETGEN_OUTPUT" ))
3982 // redirect all netgen output (mycout,myerr,cout) to _outputFileName
3983 _outputFileName = getOutputFileName();
3984 _ngcout = netgen::mycout;
3985 _ngcerr = netgen::myerr;
3986 netgen::mycout = new ofstream ( _outputFileName.c_str() );
3987 netgen::myerr = netgen::mycout;
3988 _coutBuffer = std::cout.rdbuf();
3990 cout << "NOTE: netgen output is redirected to file " << _outputFileName << endl;
3992 std::cout.rdbuf( netgen::mycout->rdbuf() );
3996 _ngMesh = Ng_NewMesh();
3999 //================================================================================
4001 * \brief Finish using netgen library
4003 //================================================================================
4005 NETGENPlugin_NetgenLibWrapper::~NETGENPlugin_NetgenLibWrapper()
4007 --instanceCounter();
4009 Ng_DeleteMesh( _ngMesh );
4013 std::cout.rdbuf( _coutBuffer );
4020 //================================================================================
4022 * \brief Set netgen mesh to delete at destruction
4024 //================================================================================
4026 void NETGENPlugin_NetgenLibWrapper::setMesh( Ng_Mesh* mesh )
4029 Ng_DeleteMesh( _ngMesh );
4033 //================================================================================
4035 * \brief Return a unique file name
4037 //================================================================================
4039 std::string NETGENPlugin_NetgenLibWrapper::getOutputFileName()
4041 std::string aTmpDir = SALOMEDS_Tool::GetTmpDir();
4043 TCollection_AsciiString aGenericName = (char*)aTmpDir.c_str();
4044 aGenericName += "NETGEN_";
4046 aGenericName += getpid();
4048 aGenericName += _getpid();
4050 aGenericName += "_";
4051 aGenericName += Abs((Standard_Integer)(long) aGenericName.ToCString());
4052 aGenericName += ".out";
4054 return aGenericName.ToCString();
4057 //================================================================================
4059 * \brief Remove "test.out" and "problemfaces" files in current directory
4061 //================================================================================
4063 void NETGENPlugin_NetgenLibWrapper::RemoveTmpFiles()
4065 bool rm = SMESH_File("test.out").remove() ;
4067 if ( rm && netgen::testout && instanceCounter() == 0 )
4069 delete netgen::testout;
4070 netgen::testout = 0;
4073 SMESH_File("problemfaces").remove();
4074 SMESH_File("occmesh.rep").remove();
4077 //================================================================================
4079 * \brief Remove file with netgen output
4081 //================================================================================
4083 void NETGENPlugin_NetgenLibWrapper::removeOutputFile()
4085 if ( !_outputFileName.empty() )
4089 delete netgen::mycout;
4090 netgen::mycout = _ngcout;
4091 netgen::myerr = _ngcerr;
4094 string tmpDir = SALOMEDS_Tool::GetDirFromPath ( _outputFileName );
4095 string aFileName = SALOMEDS_Tool::GetNameFromPath( _outputFileName ) + ".out";
4096 SALOMEDS::ListOfFileNames_var aFiles = new SALOMEDS::ListOfFileNames;
4098 aFiles[0] = aFileName.c_str();
4100 SALOMEDS_Tool::RemoveTemporaryFiles( tmpDir.c_str(), aFiles.in(), true );