1 // Copyright (C) 2007-2015 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_MeshElement.hxx>
36 #include <SMDS_MeshNode.hxx>
37 #include <SMESHDS_Mesh.hxx>
38 #include <SMESH_Block.hxx>
39 #include <SMESH_Comment.hxx>
40 #include <SMESH_ComputeError.hxx>
41 #include <SMESH_File.hxx>
42 #include <SMESH_Gen_i.hxx>
43 #include <SMESH_Mesh.hxx>
44 #include <SMESH_MesherHelper.hxx>
45 #include <SMESH_subMesh.hxx>
46 #include <StdMeshers_QuadToTriaAdaptor.hxx>
47 #include <StdMeshers_ViscousLayers2D.hxx>
49 #include <SALOMEDS_Tool.hxx>
51 #include <utilities.h>
53 #include <BRepBuilderAPI_Copy.hxx>
54 #include <BRep_Tool.hxx>
55 #include <Bnd_B3d.hxx>
56 #include <NCollection_Map.hxx>
57 #include <Standard_ErrorHandler.hxx>
58 #include <Standard_ProgramError.hxx>
59 #include <TColStd_MapOfInteger.hxx>
61 #include <TopExp_Explorer.hxx>
62 #include <TopTools_DataMapIteratorOfDataMapOfShapeInteger.hxx>
63 #include <TopTools_DataMapIteratorOfDataMapOfShapeShape.hxx>
64 #include <TopTools_DataMapOfShapeInteger.hxx>
65 #include <TopTools_DataMapOfShapeShape.hxx>
66 #include <TopTools_MapOfShape.hxx>
69 // Netgen include files
73 #include <occgeom.hpp>
74 #include <meshing.hpp>
75 //#include <ngexception.hpp>
78 extern int OCCGenerateMesh (OCCGeometry&, Mesh*&, MeshingParameters&, int, int);
80 extern int OCCGenerateMesh (OCCGeometry&, Mesh*&, int, int, char*);
82 //extern void OCCSetLocalMeshSize(OCCGeometry & geom, Mesh & mesh);
83 extern MeshingParameters mparam;
84 extern volatile multithreadt multithread;
85 extern bool merge_solids;
94 using namespace nglib;
98 #define nodeVec_ACCESS(index) ((SMDS_MeshNode*) nodeVec.at((index)))
100 #define nodeVec_ACCESS(index) ((SMDS_MeshNode*) nodeVec[index])
103 #define NGPOINT_COORDS(p) p(0),p(1),p(2)
106 // dump elements added to ng mesh
107 //#define DUMP_SEGMENTS
108 //#define DUMP_TRIANGLES
109 //#define DUMP_TRIANGLES_SCRIPT "/tmp/trias.py" //!< debug AddIntVerticesInSolids()
112 TopTools_IndexedMapOfShape ShapesWithLocalSize;
113 std::map<int,double> VertexId2LocalSize;
114 std::map<int,double> EdgeId2LocalSize;
115 std::map<int,double> FaceId2LocalSize;
117 //=============================================================================
121 //=============================================================================
123 NETGENPlugin_Mesher::NETGENPlugin_Mesher (SMESH_Mesh* mesh,
124 const TopoDS_Shape& aShape,
130 _fineness(NETGENPlugin_Hypothesis::GetDefaultFineness()),
131 _isViscousLayers2D(false),
140 SetDefaultParameters();
141 ShapesWithLocalSize.Clear();
142 VertexId2LocalSize.clear();
143 EdgeId2LocalSize.clear();
144 FaceId2LocalSize.clear();
147 //================================================================================
151 //================================================================================
153 NETGENPlugin_Mesher::~NETGENPlugin_Mesher()
161 //================================================================================
163 * Set pointer to NETGENPlugin_Mesher* field of the holder, that will be
164 * nullified at destruction of this
166 //================================================================================
168 void NETGENPlugin_Mesher::SetSelfPointer( NETGENPlugin_Mesher ** ptr )
179 //================================================================================
181 * \brief Initialize global NETGEN parameters with default values
183 //================================================================================
185 void NETGENPlugin_Mesher::SetDefaultParameters()
187 netgen::MeshingParameters& mparams = netgen::mparam;
188 // maximal mesh edge size
189 mparams.maxh = 0;//NETGENPlugin_Hypothesis::GetDefaultMaxSize();
191 // minimal number of segments per edge
192 mparams.segmentsperedge = NETGENPlugin_Hypothesis::GetDefaultNbSegPerEdge();
193 // rate of growth of size between elements
194 mparams.grading = NETGENPlugin_Hypothesis::GetDefaultGrowthRate();
195 // safety factor for curvatures (elements per radius)
196 mparams.curvaturesafety = NETGENPlugin_Hypothesis::GetDefaultNbSegPerRadius();
197 // create elements of second order
198 mparams.secondorder = NETGENPlugin_Hypothesis::GetDefaultSecondOrder();
199 // quad-dominated surface meshing
203 mparams.quad = NETGENPlugin_Hypothesis_2D::GetDefaultQuadAllowed();
204 _fineness = NETGENPlugin_Hypothesis::GetDefaultFineness();
205 mparams.uselocalh = NETGENPlugin_Hypothesis::GetDefaultSurfaceCurvature();
206 netgen::merge_solids = NETGENPlugin_Hypothesis::GetDefaultFuseEdges();
209 //=============================================================================
213 //=============================================================================
215 void SetLocalSize(TopoDS_Shape GeomShape, double LocalSize)
217 if ( GeomShape.IsNull() ) return;
218 TopAbs_ShapeEnum GeomType = GeomShape.ShapeType();
219 if (GeomType == TopAbs_COMPOUND) {
220 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()) {
221 SetLocalSize(it.Value(), LocalSize);
226 if (! ShapesWithLocalSize.Contains(GeomShape))
227 key = ShapesWithLocalSize.Add(GeomShape);
229 key = ShapesWithLocalSize.FindIndex(GeomShape);
230 if (GeomType == TopAbs_VERTEX) {
231 VertexId2LocalSize[key] = LocalSize;
232 } else if (GeomType == TopAbs_EDGE) {
233 EdgeId2LocalSize[key] = LocalSize;
234 } else if (GeomType == TopAbs_FACE) {
235 FaceId2LocalSize[key] = LocalSize;
239 //=============================================================================
241 * Pass parameters to NETGEN
243 //=============================================================================
244 void NETGENPlugin_Mesher::SetParameters(const NETGENPlugin_Hypothesis* hyp)
248 netgen::MeshingParameters& mparams = netgen::mparam;
249 // Initialize global NETGEN parameters:
250 // maximal mesh segment size
251 mparams.maxh = hyp->GetMaxSize();
252 // maximal mesh element linear size
253 mparams.minh = hyp->GetMinSize();
254 // minimal number of segments per edge
255 mparams.segmentsperedge = hyp->GetNbSegPerEdge();
256 // rate of growth of size between elements
257 mparams.grading = hyp->GetGrowthRate();
258 // safety factor for curvatures (elements per radius)
259 mparams.curvaturesafety = hyp->GetNbSegPerRadius();
260 // create elements of second order
261 mparams.secondorder = hyp->GetSecondOrder() ? 1 : 0;
262 // quad-dominated surface meshing
263 // only triangles are allowed for volumic mesh (before realizing IMP 0021676)
265 mparams.quad = hyp->GetQuadAllowed() ? 1 : 0;
266 _optimize = hyp->GetOptimize();
267 _fineness = hyp->GetFineness();
268 mparams.uselocalh = hyp->GetSurfaceCurvature();
269 netgen::merge_solids = hyp->GetFuseEdges();
272 SMESH_Gen_i* smeshGen_i = SMESH_Gen_i::GetSMESHGen();
273 CORBA::Object_var anObject = smeshGen_i->GetNS()->Resolve("/myStudyManager");
274 SALOMEDS::StudyManager_var aStudyMgr = SALOMEDS::StudyManager::_narrow(anObject);
275 SALOMEDS::Study_var myStudy = aStudyMgr->GetStudyByID(hyp->GetStudyId());
277 const NETGENPlugin_Hypothesis::TLocalSize localSizes = hyp->GetLocalSizesAndEntries();
278 NETGENPlugin_Hypothesis::TLocalSize::const_iterator it = localSizes.begin();
279 for ( ; it != localSizes.end() ; it++)
281 std::string entry = (*it).first;
282 double val = (*it).second;
284 GEOM::GEOM_Object_var aGeomObj;
285 TopoDS_Shape S = TopoDS_Shape();
286 SALOMEDS::SObject_var aSObj = myStudy->FindObjectID( entry.c_str() );
287 if (!aSObj->_is_nil()) {
288 CORBA::Object_var obj = aSObj->GetObject();
289 aGeomObj = GEOM::GEOM_Object::_narrow(obj);
292 if ( !aGeomObj->_is_nil() )
293 S = smeshGen_i->GeomObjectToShape( aGeomObj.in() );
295 SetLocalSize(S, val);
300 //=============================================================================
302 * Pass simple parameters to NETGEN
304 //=============================================================================
306 void NETGENPlugin_Mesher::SetParameters(const NETGENPlugin_SimpleHypothesis_2D* hyp)
310 SetDefaultParameters();
313 //=============================================================================
315 * Link - a pair of integer numbers
317 //=============================================================================
321 Link(int _n1, int _n2) : n1(_n1), n2(_n2) {}
322 Link() : n1(0), n2(0) {}
323 bool Contains( int n ) const { return n == n1 || n == n2; }
324 bool IsConnected( const Link& other ) const
326 return (( Contains( other.n1 ) || Contains( other.n2 )) && ( this != &other ));
330 int HashCode(const Link& aLink, int aLimit)
332 return HashCode(aLink.n1 + aLink.n2, aLimit);
335 Standard_Boolean IsEqual(const Link& aLink1, const Link& aLink2)
337 return (( aLink1.n1 == aLink2.n1 && aLink1.n2 == aLink2.n2 ) ||
338 ( aLink1.n1 == aLink2.n2 && aLink1.n2 == aLink2.n1 ));
343 //================================================================================
345 * \brief return id of netgen point corresponding to SMDS node
347 //================================================================================
348 typedef map< const SMDS_MeshNode*, int > TNode2IdMap;
350 int ngNodeId( const SMDS_MeshNode* node,
351 netgen::Mesh& ngMesh,
352 TNode2IdMap& nodeNgIdMap)
354 int newNgId = ngMesh.GetNP() + 1;
356 TNode2IdMap::iterator node_id = nodeNgIdMap.insert( make_pair( node, newNgId )).first;
358 if ( node_id->second == newNgId)
360 #if defined(DUMP_SEGMENTS) || defined(DUMP_TRIANGLES)
361 cout << "Ng " << newNgId << " - " << node;
363 netgen::MeshPoint p( netgen::Point<3> (node->X(), node->Y(), node->Z()) );
364 ngMesh.AddPoint( p );
366 return node_id->second;
369 //================================================================================
371 * \brief Return computed EDGEs connected to the given one
373 //================================================================================
375 list< TopoDS_Edge > getConnectedEdges( const TopoDS_Edge& edge,
376 const TopoDS_Face& face,
377 const set< SMESH_subMesh* > & computedSM,
378 const SMESH_MesherHelper& helper,
379 map< SMESH_subMesh*, set< int > >& addedEdgeSM2Faces)
382 list< TopoDS_Edge > edges;
383 list< int > nbEdgesInWire;
384 /*int nbWires =*/ SMESH_Block::GetOrderedEdges( face, edges, nbEdgesInWire);
386 // find <edge> within <edges>
387 list< TopoDS_Edge >::iterator eItFwd = edges.begin();
388 for ( ; eItFwd != edges.end(); ++eItFwd )
389 if ( edge.IsSame( *eItFwd ))
391 if ( eItFwd == edges.end()) return list< TopoDS_Edge>();
393 if ( eItFwd->Orientation() >= TopAbs_INTERNAL )
395 // connected INTERNAL edges returned from GetOrderedEdges() are wrongly oriented
396 // so treat each INTERNAL edge separately
397 TopoDS_Edge e = *eItFwd;
399 edges.push_back( e );
403 // get all computed EDGEs connected to <edge>
405 list< TopoDS_Edge >::iterator eItBack = eItFwd, ePrev;
406 TopoDS_Vertex vCommon;
407 TopTools_MapOfShape eAdded; // map used not to add a seam edge twice to <edges>
410 // put edges before <edge> to <edges> back
411 while ( edges.begin() != eItFwd )
412 edges.splice( edges.end(), edges, edges.begin() );
416 while ( ++eItFwd != edges.end() )
418 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( *eItFwd );
420 bool connected = TopExp::CommonVertex( *ePrev, *eItFwd, vCommon );
421 bool computed = sm->IsMeshComputed();
422 bool added = addedEdgeSM2Faces[ sm ].count( helper.GetSubShapeID() );
423 bool doubled = !eAdded.Add( *eItFwd );
424 bool orientOK = (( ePrev ->Orientation() < TopAbs_INTERNAL ) ==
425 ( eItFwd->Orientation() < TopAbs_INTERNAL ) );
426 if ( !connected || !computed || !orientOK || added || doubled )
428 // stop advancement; move edges from tail to head
429 while ( edges.back() != *ePrev )
430 edges.splice( edges.begin(), edges, --edges.end() );
436 while ( eItBack != edges.begin() )
440 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( *eItBack );
442 bool connected = TopExp::CommonVertex( *ePrev, *eItBack, vCommon );
443 bool computed = sm->IsMeshComputed();
444 bool added = addedEdgeSM2Faces[ sm ].count( helper.GetSubShapeID() );
445 bool doubled = !eAdded.Add( *eItBack );
446 bool orientOK = (( ePrev ->Orientation() < TopAbs_INTERNAL ) ==
447 ( eItBack->Orientation() < TopAbs_INTERNAL ) );
448 if ( !connected || !computed || !orientOK || added || doubled)
451 edges.erase( edges.begin(), ePrev );
455 if ( edges.front() != edges.back() )
457 // assure that the 1st vertex is meshed
458 TopoDS_Edge eLast = edges.back();
459 while ( !SMESH_Algo::VertexNode( SMESH_MesherHelper::IthVertex( 0, edges.front()), helper.GetMeshDS())
461 edges.front() != eLast )
462 edges.splice( edges.end(), edges, edges.begin() );
467 //================================================================================
469 * \brief Make triangulation of a shape precise enough
471 //================================================================================
473 void updateTriangulation( const TopoDS_Shape& shape )
475 // static set< Poly_Triangulation* > updated;
477 // TopLoc_Location loc;
478 // TopExp_Explorer fExp( shape, TopAbs_FACE );
479 // for ( ; fExp.More(); fExp.Next() )
481 // Handle(Poly_Triangulation) triangulation =
482 // BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
483 // if ( triangulation.IsNull() ||
484 // updated.insert( triangulation.operator->() ).second )
486 // BRepTools::Clean (shape);
489 BRepMesh_IncrementalMesh e(shape, 0.01, true);
491 catch (Standard_Failure)
494 // updated.erase( triangulation.operator->() );
495 // triangulation = BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
496 // updated.insert( triangulation.operator->() );
500 //================================================================================
502 * \brief Returns a medium node either existing in SMESH of created by NETGEN
503 * \param [in] corner1 - corner node 1
504 * \param [in] corner2 - corner node 2
505 * \param [in] defaultMedium - the node created by NETGEN
506 * \param [in] helper - holder of medium nodes existing in SMESH
507 * \return const SMDS_MeshNode* - the result node
509 //================================================================================
511 const SMDS_MeshNode* mediumNode( const SMDS_MeshNode* corner1,
512 const SMDS_MeshNode* corner2,
513 const SMDS_MeshNode* defaultMedium,
514 const SMESH_MesherHelper* helper)
518 TLinkNodeMap::const_iterator l2n =
519 helper->GetTLinkNodeMap().find( SMESH_TLink( corner1, corner2 ));
520 if ( l2n != helper->GetTLinkNodeMap().end() )
521 defaultMedium = l2n->second;
523 return defaultMedium;
526 //================================================================================
528 * \brief Assure that mesh on given shapes is quadratic
530 //================================================================================
532 void makeQuadratic( const TopTools_IndexedMapOfShape& shapes,
535 for ( int i = 1; i <= shapes.Extent(); ++i )
537 SMESHDS_SubMesh* smDS = mesh->GetMeshDS()->MeshElements( shapes(i) );
538 if ( !smDS ) continue;
539 SMDS_ElemIteratorPtr elemIt = smDS->GetElements();
540 if ( !elemIt->more() ) continue;
541 const SMDS_MeshElement* e = elemIt->next();
542 if ( !e || e->IsQuadratic() )
545 TIDSortedElemSet elems;
547 while ( elemIt->more() )
548 elems.insert( elems.end(), elemIt->next() );
550 SMESH_MeshEditor( mesh ).ConvertToQuadratic( /*3d=*/false, elems, /*biQuad=*/false );
556 //================================================================================
558 * \brief Initialize netgen::OCCGeometry with OCCT shape
560 //================================================================================
562 void NETGENPlugin_Mesher::PrepareOCCgeometry(netgen::OCCGeometry& occgeo,
563 const TopoDS_Shape& shape,
565 list< SMESH_subMesh* > * meshedSM,
566 NETGENPlugin_Internals* intern)
568 updateTriangulation( shape );
571 BRepBndLib::Add (shape, bb);
572 double x1,y1,z1,x2,y2,z2;
573 bb.Get (x1,y1,z1,x2,y2,z2);
574 MESSAGE("shape bounding box:\n" <<
575 "(" << x1 << " " << y1 << " " << z1 << ") " <<
576 "(" << x2 << " " << y2 << " " << z2 << ")");
577 netgen::Point<3> p1 = netgen::Point<3> (x1,y1,z1);
578 netgen::Point<3> p2 = netgen::Point<3> (x2,y2,z2);
579 occgeo.boundingbox = netgen::Box<3> (p1,p2);
581 occgeo.shape = shape;
584 // fill maps of shapes of occgeo with not yet meshed subshapes
586 // get root submeshes
587 list< SMESH_subMesh* > rootSM;
588 const int shapeID = mesh.GetMeshDS()->ShapeToIndex( shape );
589 if ( shapeID > 0 ) { // SMESH_subMesh with ID 0 may exist, don't use it!
590 rootSM.push_back( mesh.GetSubMesh( shape ));
593 for ( TopoDS_Iterator it( shape ); it.More(); it.Next() )
594 rootSM.push_back( mesh.GetSubMesh( it.Value() ));
597 // add subshapes of empty submeshes
598 list< SMESH_subMesh* >::iterator rootIt = rootSM.begin(), rootEnd = rootSM.end();
599 for ( ; rootIt != rootEnd; ++rootIt ) {
600 SMESH_subMesh * root = *rootIt;
601 SMESH_subMeshIteratorPtr smIt = root->getDependsOnIterator(/*includeSelf=*/true,
602 /*complexShapeFirst=*/true);
603 // to find a right orientation of subshapes (PAL20462)
604 TopTools_IndexedMapOfShape subShapes;
605 TopExp::MapShapes(root->GetSubShape(), subShapes);
606 while ( smIt->more() )
608 SMESH_subMesh* sm = smIt->next();
609 TopoDS_Shape shape = sm->GetSubShape();
610 if ( intern && intern->isShapeToPrecompute( shape ))
612 if ( !meshedSM || sm->IsEmpty() )
614 if ( shape.ShapeType() != TopAbs_VERTEX )
615 shape = subShapes( subShapes.FindIndex( shape ));// shape -> index -> oriented shape
616 if ( shape.Orientation() >= TopAbs_INTERNAL )
617 shape.Orientation( TopAbs_FORWARD ); // isuue 0020676
618 switch ( shape.ShapeType() ) {
619 case TopAbs_FACE : occgeo.fmap.Add( shape ); break;
620 case TopAbs_EDGE : occgeo.emap.Add( shape ); break;
621 case TopAbs_VERTEX: occgeo.vmap.Add( shape ); break;
622 case TopAbs_SOLID :occgeo.somap.Add( shape ); break;
626 // collect submeshes of meshed shapes
629 const int dim = SMESH_Gen::GetShapeDim( shape );
630 meshedSM[ dim ].push_back( sm );
634 occgeo.facemeshstatus.SetSize (occgeo.fmap.Extent());
635 occgeo.facemeshstatus = 0;
636 occgeo.face_maxh_modified.SetSize(occgeo.fmap.Extent());
637 occgeo.face_maxh_modified = 0;
638 occgeo.face_maxh.SetSize(occgeo.fmap.Extent());
639 occgeo.face_maxh = netgen::mparam.maxh;
642 //================================================================================
644 * \brief Return a default min size value suitable for the given geometry.
646 //================================================================================
648 double NETGENPlugin_Mesher::GetDefaultMinSize(const TopoDS_Shape& geom,
649 const double maxSize)
651 updateTriangulation( geom );
655 const int* pi[4] = { &i1, &i2, &i3, &i1 };
658 TopExp_Explorer fExp( geom, TopAbs_FACE );
659 for ( ; fExp.More(); fExp.Next() )
661 Handle(Poly_Triangulation) triangulation =
662 BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
663 if ( triangulation.IsNull() ) continue;
664 const double fTol = BRep_Tool::Tolerance( TopoDS::Face( fExp.Current() ));
665 const TColgp_Array1OfPnt& points = triangulation->Nodes();
666 const Poly_Array1OfTriangle& trias = triangulation->Triangles();
667 for ( int iT = trias.Lower(); iT <= trias.Upper(); ++iT )
669 trias(iT).Get( i1, i2, i3 );
670 for ( int j = 0; j < 3; ++j )
672 double dist2 = points(*pi[j]).SquareDistance( points( *pi[j+1] ));
673 if ( dist2 < minh && fTol*fTol < dist2 )
675 bb.Add( points(*pi[j]));
679 if ( minh > 0.25 * bb.SquareExtent() ) // simple geometry, rough triangulation
681 minh = 1e-3 * sqrt( bb.SquareExtent());
682 //cout << "BND BOX minh = " <<minh << endl;
686 minh = 3 * sqrt( minh ); // triangulation for visualization is rather fine
687 //cout << "TRIANGULATION minh = " <<minh << endl;
689 if ( minh > 0.5 * maxSize )
695 //================================================================================
697 * \brief Restrict size of elements at a given point
699 //================================================================================
701 void NETGENPlugin_Mesher::RestrictLocalSize(netgen::Mesh& ngMesh,
704 const bool overrideMinH)
706 if ( size <= std::numeric_limits<double>::min() )
708 if ( netgen::mparam.minh > size )
712 ngMesh.SetMinimalH( size );
713 netgen::mparam.minh = size;
717 size = netgen::mparam.minh;
720 netgen::Point3d pi(p.X(), p.Y(), p.Z());
721 ngMesh.RestrictLocalH( pi, size );
724 //================================================================================
726 * \brief fill ngMesh with nodes and elements of computed submeshes
728 //================================================================================
730 bool NETGENPlugin_Mesher::FillNgMesh(netgen::OCCGeometry& occgeom,
731 netgen::Mesh& ngMesh,
732 vector<const SMDS_MeshNode*>& nodeVec,
733 const list< SMESH_subMesh* > & meshedSM,
734 SMESH_MesherHelper* quadHelper,
735 SMESH_ProxyMesh::Ptr proxyMesh)
737 TNode2IdMap nodeNgIdMap;
738 for ( size_t i = 1; i < nodeVec.size(); ++i )
739 nodeNgIdMap.insert( make_pair( nodeVec[i], i ));
741 TopTools_MapOfShape visitedShapes;
742 map< SMESH_subMesh*, set< int > > visitedEdgeSM2Faces;
743 set< SMESH_subMesh* > computedSM( meshedSM.begin(), meshedSM.end() );
745 SMESH_MesherHelper helper (*_mesh);
747 int faceNgID = ngMesh.GetNFD();
749 list< SMESH_subMesh* >::const_iterator smIt, smEnd = meshedSM.end();
750 for ( smIt = meshedSM.begin(); smIt != smEnd; ++smIt )
752 SMESH_subMesh* sm = *smIt;
753 if ( !visitedShapes.Add( sm->GetSubShape() ))
756 const SMESHDS_SubMesh * smDS = sm->GetSubMeshDS();
757 if ( !smDS ) continue;
759 switch ( sm->GetSubShape().ShapeType() )
761 case TopAbs_EDGE: { // EDGE
762 // ----------------------
763 TopoDS_Edge geomEdge = TopoDS::Edge( sm->GetSubShape() );
764 if ( geomEdge.Orientation() >= TopAbs_INTERNAL )
765 geomEdge.Orientation( TopAbs_FORWARD ); // issue 0020676
767 // Add ng segments for each not meshed FACE the EDGE bounds
768 PShapeIteratorPtr fIt = helper.GetAncestors( geomEdge, *sm->GetFather(), TopAbs_FACE );
769 while ( const TopoDS_Shape * anc = fIt->next() )
771 faceNgID = occgeom.fmap.FindIndex( *anc );
773 continue; // meshed face
775 int faceSMDSId = helper.GetMeshDS()->ShapeToIndex( *anc );
776 if ( visitedEdgeSM2Faces[ sm ].count( faceSMDSId ))
777 continue; // already treated EDGE
779 TopoDS_Face face = TopoDS::Face( occgeom.fmap( faceNgID ));
780 if ( face.Orientation() >= TopAbs_INTERNAL )
781 face.Orientation( TopAbs_FORWARD ); // issue 0020676
783 // get all meshed EDGEs of the FACE connected to geomEdge (issue 0021140)
784 helper.SetSubShape( face );
785 list< TopoDS_Edge > edges = getConnectedEdges( geomEdge, face, computedSM, helper,
786 visitedEdgeSM2Faces );
788 continue; // wrong ancestor?
790 // find out orientation of <edges> within <face>
791 TopoDS_Edge eNotSeam = edges.front();
792 if ( helper.HasSeam() )
794 list< TopoDS_Edge >::iterator eIt = edges.begin();
795 while ( helper.IsRealSeam( *eIt )) ++eIt;
796 if ( eIt != edges.end() )
799 TopAbs_Orientation fOri = helper.GetSubShapeOri( face, eNotSeam );
800 bool isForwad = ( fOri == eNotSeam.Orientation() || fOri >= TopAbs_INTERNAL );
802 // get all nodes from connected <edges>
803 const bool isQuad = smDS->IsQuadratic();
804 StdMeshers_FaceSide fSide( face, edges, _mesh, isForwad, isQuad );
805 const vector<UVPtStruct>& points = fSide.GetUVPtStruct();
806 if ( points.empty() )
807 return false; // invalid node params?
808 int i, nbSeg = fSide.NbSegments();
810 // remember EDGEs of fSide to treat only once
811 for ( int iE = 0; iE < fSide.NbEdges(); ++iE )
812 visitedEdgeSM2Faces[ helper.GetMesh()->GetSubMesh( fSide.Edge(iE )) ].insert(faceSMDSId);
814 double otherSeamParam = 0;
819 int prevNgId = ngNodeId( points[0].node, ngMesh, nodeNgIdMap );
821 for ( i = 0; i < nbSeg; ++i )
823 const UVPtStruct& p1 = points[ i ];
824 const UVPtStruct& p2 = points[ i+1 ];
826 if ( p1.node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX ) //an EDGE begins
829 if ( helper.IsRealSeam( p1.node->getshapeId() ))
831 TopoDS_Edge e = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
832 isSeam = helper.IsRealSeam( e );
835 otherSeamParam = helper.GetOtherParam( helper.GetPeriodicIndex() & 1 ? p2.u : p2.v );
842 seg[1] = prevNgId = ngNodeId( p2.node, ngMesh, nodeNgIdMap );
843 // node param on curve
844 seg.epgeominfo[ 0 ].dist = p1.param;
845 seg.epgeominfo[ 1 ].dist = p2.param;
847 seg.epgeominfo[ 0 ].u = p1.u;
848 seg.epgeominfo[ 0 ].v = p1.v;
849 seg.epgeominfo[ 1 ].u = p2.u;
850 seg.epgeominfo[ 1 ].v = p2.v;
852 //geomEdge = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
853 //seg.epgeominfo[ 0 ].edgenr = seg.epgeominfo[ 1 ].edgenr = occgeom.emap.FindIndex( geomEdge );
855 //seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
856 seg.si = faceNgID; // = geom.fmap.FindIndex (face);
857 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
858 ngMesh.AddSegment (seg);
860 SMESH_TNodeXYZ np1( p1.node ), np2( p2.node );
861 RestrictLocalSize( ngMesh, 0.5*(np1+np2), (np1-np2).Modulus() );
864 cout << "Segment: " << seg.edgenr << " on SMESH face " << helper.GetMeshDS()->ShapeToIndex( face ) << endl
865 << "\tface index: " << seg.si << endl
866 << "\tp1: " << seg[0] << endl
867 << "\tp2: " << seg[1] << endl
868 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
869 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
870 //<< "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
871 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
872 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl;
873 //<< "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
877 if ( helper.GetPeriodicIndex() && 1 ) {
878 seg.epgeominfo[ 0 ].u = otherSeamParam;
879 seg.epgeominfo[ 1 ].u = otherSeamParam;
880 swap (seg.epgeominfo[0].v, seg.epgeominfo[1].v);
882 seg.epgeominfo[ 0 ].v = otherSeamParam;
883 seg.epgeominfo[ 1 ].v = otherSeamParam;
884 swap (seg.epgeominfo[0].u, seg.epgeominfo[1].u);
886 swap (seg[0], seg[1]);
887 swap (seg.epgeominfo[0].dist, seg.epgeominfo[1].dist);
888 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
889 ngMesh.AddSegment (seg);
891 cout << "Segment: " << seg.edgenr << endl
892 << "\t is SEAM (reverse) of the previous. "
893 << " Other " << (helper.GetPeriodicIndex() && 1 ? "U" : "V")
894 << " = " << otherSeamParam << endl;
897 else if ( fOri == TopAbs_INTERNAL )
899 swap (seg[0], seg[1]);
900 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
901 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
902 ngMesh.AddSegment (seg);
904 cout << "Segment: " << seg.edgenr << endl << "\t is REVERSE of the previous" << endl;
908 } // loop on geomEdge ancestors
910 if ( quadHelper ) // remember medium nodes of sub-meshes
912 SMDS_ElemIteratorPtr edges = smDS->GetElements();
913 while ( edges->more() )
915 const SMDS_MeshElement* e = edges->next();
916 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshEdge*>( e )))
922 } // case TopAbs_EDGE
924 case TopAbs_FACE: { // FACE
925 // ----------------------
926 const TopoDS_Face& geomFace = TopoDS::Face( sm->GetSubShape() );
927 helper.SetSubShape( geomFace );
928 bool isInternalFace = ( geomFace.Orientation() == TopAbs_INTERNAL );
930 // Find solids the geomFace bounds
931 int solidID1 = 0, solidID2 = 0;
932 StdMeshers_QuadToTriaAdaptor* quadAdaptor =
933 dynamic_cast<StdMeshers_QuadToTriaAdaptor*>( proxyMesh.get() );
936 solidID1 = occgeom.somap.FindIndex( quadAdaptor->GetShape() );
940 PShapeIteratorPtr solidIt = helper.GetAncestors( geomFace, *sm->GetFather(), TopAbs_SOLID);
941 while ( const TopoDS_Shape * solid = solidIt->next() )
943 int id = occgeom.somap.FindIndex ( *solid );
944 if ( solidID1 && id != solidID1 ) solidID2 = id;
948 // Add ng face descriptors of meshed faces
950 ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(faceNgID, solidID1, solidID2, 0));
952 // if second oreder is required, even already meshed faces must be passed to NETGEN
953 int fID = occgeom.fmap.Add( geomFace );
954 while ( fID < faceNgID ) // geomFace is already in occgeom.fmap, add a copy
955 fID = occgeom.fmap.Add( BRepBuilderAPI_Copy( geomFace, /*copyGeom=*/false ));
956 // Problem with the second order in a quadrangular mesh remains.
957 // 1) All quadrangles generated by NETGEN are moved to an inexistent face
958 // by FillSMesh() (find "AddFaceDescriptor")
959 // 2) Temporary triangles generated by StdMeshers_QuadToTriaAdaptor
960 // are on faces where quadrangles were.
961 // Due to these 2 points, wrong geom faces are used while conversion to qudratic
962 // of the mentioned above quadrangles and triangles
964 // Orient the face correctly in solidID1 (issue 0020206)
965 bool reverse = false;
967 TopoDS_Shape solid = occgeom.somap( solidID1 );
968 TopAbs_Orientation faceOriInSolid = helper.GetSubShapeOri( solid, geomFace );
969 if ( faceOriInSolid >= 0 )
971 helper.IsReversedSubMesh( TopoDS::Face( geomFace.Oriented( faceOriInSolid )));
974 // Add surface elements
976 netgen::Element2d tri(3);
977 tri.SetIndex ( faceNgID );
978 SMESH_TNodeXYZ xyz[3];
980 #ifdef DUMP_TRIANGLES
981 cout << "SMESH face " << helper.GetMeshDS()->ShapeToIndex( geomFace )
982 << " internal="<<isInternalFace << endl;
985 smDS = proxyMesh->GetSubMesh( geomFace );
987 SMDS_ElemIteratorPtr faces = smDS->GetElements();
988 while ( faces->more() )
990 const SMDS_MeshElement* f = faces->next();
991 if ( f->NbNodes() % 3 != 0 ) // not triangle
993 PShapeIteratorPtr solidIt=helper.GetAncestors(geomFace,*sm->GetFather(),TopAbs_SOLID);
994 if ( const TopoDS_Shape * solid = solidIt->next() )
995 sm = _mesh->GetSubMesh( *solid );
996 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
997 smError.reset( new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH,"Not triangle sub-mesh"));
998 smError->myBadElements.push_back( f );
1002 for ( int i = 0; i < 3; ++i )
1004 const SMDS_MeshNode* node = f->GetNode( i ), * inFaceNode=0;
1007 // get node UV on face
1008 int shapeID = node->getshapeId();
1009 if ( helper.IsSeamShape( shapeID ))
1011 if ( helper.IsSeamShape( f->GetNodeWrap( i+1 )->getshapeId() ))
1012 inFaceNode = f->GetNodeWrap( i-1 );
1014 inFaceNode = f->GetNodeWrap( i+1 );
1016 gp_XY uv = helper.GetNodeUV( geomFace, node, inFaceNode );
1018 int ind = reverse ? 3-i : i+1;
1019 tri.GeomInfoPi(ind).u = uv.X();
1020 tri.GeomInfoPi(ind).v = uv.Y();
1021 tri.PNum (ind) = ngNodeId( node, ngMesh, nodeNgIdMap );
1024 // pass a triangle size to NG size-map
1025 double size = ( ( xyz[0] - xyz[1] ).Modulus() +
1026 ( xyz[1] - xyz[2] ).Modulus() +
1027 ( xyz[2] - xyz[0] ).Modulus() ) / 3;
1028 gp_XYZ gc = ( xyz[0] + xyz[1] + xyz[2] ) / 3;
1029 RestrictLocalSize( ngMesh, gc, size, /*overrideMinH=*/false );
1031 ngMesh.AddSurfaceElement (tri);
1032 #ifdef DUMP_TRIANGLES
1033 cout << tri << endl;
1036 if ( isInternalFace )
1038 swap( tri[1], tri[2] );
1039 ngMesh.AddSurfaceElement (tri);
1040 #ifdef DUMP_TRIANGLES
1041 cout << tri << endl;
1046 if ( quadHelper ) // remember medium nodes of sub-meshes
1048 SMDS_ElemIteratorPtr faces = smDS->GetElements();
1049 while ( faces->more() )
1051 const SMDS_MeshElement* f = faces->next();
1052 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshFace*>( f )))
1058 } // case TopAbs_FACE
1060 case TopAbs_VERTEX: { // VERTEX
1061 // --------------------------
1062 // issue 0021405. Add node only if a VERTEX is shared by a not meshed EDGE,
1063 // else netgen removes a free node and nodeVector becomes invalid
1064 PShapeIteratorPtr ansIt = helper.GetAncestors( sm->GetSubShape(),
1068 while ( const TopoDS_Shape* e = ansIt->next() )
1070 SMESH_subMesh* eSub = helper.GetMesh()->GetSubMesh( *e );
1071 if (( toAdd = ( eSub->IsEmpty() && !SMESH_Algo::isDegenerated( TopoDS::Edge( *e )))))
1076 SMDS_NodeIteratorPtr nodeIt = smDS->GetNodes();
1077 if ( nodeIt->more() )
1078 ngNodeId( nodeIt->next(), ngMesh, nodeNgIdMap );
1084 } // loop on submeshes
1087 nodeVec.resize( ngMesh.GetNP() + 1 );
1088 TNode2IdMap::iterator node_NgId, nodeNgIdEnd = nodeNgIdMap.end();
1089 for ( node_NgId = nodeNgIdMap.begin(); node_NgId != nodeNgIdEnd; ++node_NgId)
1090 nodeVec[ node_NgId->second ] = node_NgId->first;
1095 //================================================================================
1097 * \brief Duplicate mesh faces on internal geom faces
1099 //================================================================================
1101 void NETGENPlugin_Mesher::FixIntFaces(const netgen::OCCGeometry& occgeom,
1102 netgen::Mesh& ngMesh,
1103 NETGENPlugin_Internals& internalShapes)
1105 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1107 // find ng indices of internal faces
1109 for ( int ngFaceID = 1; ngFaceID <= occgeom.fmap.Extent(); ++ngFaceID )
1111 int smeshID = meshDS->ShapeToIndex( occgeom.fmap( ngFaceID ));
1112 if ( internalShapes.isInternalShape( smeshID ))
1113 ngFaceIds.insert( ngFaceID );
1115 if ( !ngFaceIds.empty() )
1118 int i, nbFaces = ngMesh.GetNSE();
1119 for ( i = 1; i <= nbFaces; ++i)
1121 netgen::Element2d elem = ngMesh.SurfaceElement(i);
1122 if ( ngFaceIds.count( elem.GetIndex() ))
1124 swap( elem[1], elem[2] );
1125 ngMesh.AddSurfaceElement (elem);
1131 //================================================================================
1133 * \brief Tries to heal the mesh on a FACE. The FACE is supposed to be partially
1134 * meshed due to NETGEN failure
1135 * \param [in] occgeom - geometry
1136 * \param [in,out] ngMesh - the mesh to fix
1137 * \param [inout] faceID - ID of the FACE to fix the mesh on
1138 * \return bool - is mesh is or becomes OK
1140 //================================================================================
1142 bool NETGENPlugin_Mesher::FixFaceMesh(const netgen::OCCGeometry& occgeom,
1143 netgen::Mesh& ngMesh,
1146 // we address a case where the FACE is almost fully meshed except small holes
1147 // of usually triangular shape at FACE boundary (IPAL52861)
1149 // The case appeared to be not simple: holes only look triangular but
1150 // indeed are a self intersecting polygon. A reason of the bug was in coincident
1151 // NG points on a seam edge. But the code below is very nice, leave it for
1156 if ( occgeom.fmap.Extent() < faceID )
1158 //const TopoDS_Face& face = TopoDS::Face( occgeom.fmap( faceID ));
1160 // find free links on the FACE
1161 NCollection_Map<Link> linkMap;
1162 for ( int iF = 1; iF <= ngMesh.GetNSE(); ++iF )
1164 const netgen::Element2d& elem = ngMesh.SurfaceElement(iF);
1165 if ( faceID != elem.GetIndex() )
1167 int n0 = elem[ elem.GetNP() - 1 ];
1168 for ( int i = 0; i < elem.GetNP(); ++i )
1171 Link link( n0, n1 );
1172 if ( !linkMap.Add( link ))
1173 linkMap.Remove( link );
1177 // add/remove boundary links
1178 for ( int iSeg = 1; iSeg <= ngMesh.GetNSeg(); ++iSeg )
1180 const netgen::Segment& seg = ngMesh.LineSegment( iSeg );
1181 if ( seg.si != faceID ) // !edgeIDs.Contains( seg.edgenr ))
1183 Link link( seg[1], seg[0] ); // reverse!!!
1184 if ( !linkMap.Add( link ))
1185 linkMap.Remove( link );
1187 if ( linkMap.IsEmpty() )
1189 if ( linkMap.Extent() < 3 )
1192 // make triangles of the links
1194 netgen::Element2d tri(3);
1195 tri.SetIndex ( faceID );
1197 NCollection_Map<Link>::Iterator linkIt( linkMap );
1198 Link link1 = linkIt.Value();
1199 // look for a link connected to link1
1200 NCollection_Map<Link>::Iterator linkIt2 = linkIt;
1201 for ( linkIt2.Next(); linkIt2.More(); linkIt2.Next() )
1203 const Link& link2 = linkIt2.Value();
1204 if ( link2.IsConnected( link1 ))
1206 // look for a link connected to both link1 and link2
1207 NCollection_Map<Link>::Iterator linkIt3 = linkIt2;
1208 for ( linkIt3.Next(); linkIt3.More(); linkIt3.Next() )
1210 const Link& link3 = linkIt3.Value();
1211 if ( link3.IsConnected( link1 ) &&
1212 link3.IsConnected( link2 ) )
1217 tri[2] = ( link2.Contains( link1.n1 ) ? link2.n1 : link3.n1 );
1218 if ( tri[0] == tri[2] || tri[1] == tri[2] )
1220 ngMesh.AddSurfaceElement( tri );
1222 // prepare for the next tria search
1223 if ( linkMap.Extent() == 3 )
1225 linkMap.Remove( link3 );
1226 linkMap.Remove( link2 );
1228 linkMap.Remove( link1 );
1229 link1 = linkIt.Value();
1242 //================================================================================
1243 // define gp_XY_Subtracted pointer to function calling gp_XY::Subtracted(gp_XY)
1244 gp_XY_FunPtr(Subtracted);
1245 //gp_XY_FunPtr(Added);
1247 //================================================================================
1249 * \brief Evaluate distance between two 2d points along the surface
1251 //================================================================================
1253 double evalDist( const gp_XY& uv1,
1255 const Handle(Geom_Surface)& surf,
1256 const int stopHandler=-1)
1258 if ( stopHandler > 0 ) // continue recursion
1260 gp_XY mid = SMESH_MesherHelper::GetMiddleUV( surf, uv1, uv2 );
1261 return evalDist( uv1,mid, surf, stopHandler-1 ) + evalDist( mid,uv2, surf, stopHandler-1 );
1263 double dist3D = surf->Value( uv1.X(), uv1.Y() ).Distance( surf->Value( uv2.X(), uv2.Y() ));
1264 if ( stopHandler == 0 ) // stop recursion
1267 // start recursion if necessary
1268 double dist2D = SMESH_MesherHelper::ApplyIn2D(surf, uv1, uv2, gp_XY_Subtracted, 0).Modulus();
1269 if ( fabs( dist3D - dist2D ) < dist2D * 1e-10 )
1270 return dist3D; // equal parametrization of a planar surface
1272 return evalDist( uv1, uv2, surf, 3 ); // start recursion
1275 //================================================================================
1277 * \brief Data of vertex internal in geom face
1279 //================================================================================
1283 gp_XY uv; //!< UV in face parametric space
1284 int ngId; //!< ng id of corrsponding node
1285 gp_XY uvClose; //!< UV of closest boundary node
1286 int ngIdClose; //!< ng id of closest boundary node
1289 //================================================================================
1291 * \brief Data of vertex internal in solid
1293 //================================================================================
1297 int ngId; //!< ng id of corresponding node
1298 int ngIdClose; //!< ng id of closest 2d mesh element
1299 int ngIdCloseN; //!< ng id of closest node of the closest 2d mesh element
1302 inline double dist2(const netgen::MeshPoint& p1, const netgen::MeshPoint& p2)
1304 return gp_Pnt( NGPOINT_COORDS(p1)).SquareDistance( gp_Pnt( NGPOINT_COORDS(p2)));
1308 //================================================================================
1310 * \brief Make netgen take internal vertices in faces into account by adding
1311 * segments including internal vertices
1313 * This function works in supposition that 1D mesh is already computed in ngMesh
1315 //================================================================================
1317 void NETGENPlugin_Mesher::AddIntVerticesInFaces(const netgen::OCCGeometry& occgeom,
1318 netgen::Mesh& ngMesh,
1319 vector<const SMDS_MeshNode*>& nodeVec,
1320 NETGENPlugin_Internals& internalShapes)
1322 if ((int) nodeVec.size() < ngMesh.GetNP() )
1323 nodeVec.resize( ngMesh.GetNP(), 0 );
1325 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1326 SMESH_MesherHelper helper( internalShapes.getMesh() );
1328 const map<int,list<int> >& face2Vert = internalShapes.getFacesWithVertices();
1329 map<int,list<int> >::const_iterator f2v = face2Vert.begin();
1330 for ( ; f2v != face2Vert.end(); ++f2v )
1332 const TopoDS_Face& face = TopoDS::Face( meshDS->IndexToShape( f2v->first ));
1333 if ( face.IsNull() ) continue;
1334 int faceNgID = occgeom.fmap.FindIndex (face);
1335 if ( faceNgID < 0 ) continue;
1337 TopLoc_Location loc;
1338 Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc);
1340 helper.SetSubShape( face );
1341 helper.SetElementsOnShape( true );
1343 // Get data of internal vertices and add them to ngMesh
1345 multimap< double, TIntVData > dist2VData; // sort vertices by distance from boundary nodes
1347 int i, nbSegInit = ngMesh.GetNSeg();
1349 // boundary characteristics
1350 double totSegLen2D = 0;
1353 const list<int>& iVertices = f2v->second;
1354 list<int>::const_iterator iv = iVertices.begin();
1355 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1358 // get node on vertex
1359 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1360 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1363 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1364 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1365 nV = SMESH_Algo::VertexNode( V, meshDS );
1366 if ( !nV ) continue;
1369 netgen::MeshPoint mp( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1370 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1371 vData.ngId = ngMesh.GetNP();
1372 nodeVec.push_back( nV );
1376 vData.uv = helper.GetNodeUV( face, nV, 0, &uvOK );
1377 if ( !uvOK ) helper.CheckNodeUV( face, nV, vData.uv, BRep_Tool::Tolerance(V),/*force=*/1);
1379 // loop on all segments of the face to find the node closest to vertex and to count
1380 // average segment 2d length
1381 double closeDist2 = numeric_limits<double>::max(), dist2;
1383 for (i = 1; i <= ngMesh.GetNSeg(); ++i)
1385 netgen::Segment & seg = ngMesh.LineSegment(i);
1386 if ( seg.si != faceNgID ) continue;
1388 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1390 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1391 if ( ngIdLast == seg[ iEnd ] ) continue;
1392 dist2 = helper.ApplyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1393 if ( dist2 < closeDist2 )
1394 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1395 ngIdLast = seg[ iEnd ];
1399 totSegLen2D += helper.ApplyIn2D(surf, uv[0], uv[1], gp_XY_Subtracted, false).Modulus();
1403 dist2VData.insert( make_pair( closeDist2, vData ));
1406 if ( totNbSeg == 0 ) break;
1407 double avgSegLen2d = totSegLen2D / totNbSeg;
1409 // Loop on vertices to add segments
1411 multimap< double, TIntVData >::iterator dist_vData = dist2VData.begin();
1412 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1414 double closeDist2 = dist_vData->first, dist2;
1415 TIntVData & vData = dist_vData->second;
1417 // try to find more close node among segments added for internal vertices
1418 for (i = nbSegInit+1; i <= ngMesh.GetNSeg(); ++i)
1420 netgen::Segment & seg = ngMesh.LineSegment(i);
1421 if ( seg.si != faceNgID ) continue;
1423 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1425 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1426 dist2 = helper.ApplyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1427 if ( dist2 < closeDist2 )
1428 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1431 // decide whether to use the closest node as the second end of segment or to
1432 // create a new point
1433 int segEnd1 = vData.ngId;
1434 int segEnd2 = vData.ngIdClose; // to use closest node
1435 gp_XY uvV = vData.uv, uvP = vData.uvClose;
1436 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1437 double nodeDist2D = sqrt( closeDist2 );
1438 double nodeDist3D = evalDist( vData.uv, vData.uvClose, surf );
1439 bool avgLenOK = ( avgSegLen2d < 0.75 * nodeDist2D );
1440 bool hintLenOK = ( segLenHint < 0.75 * nodeDist3D );
1441 //cout << "uvV " << uvV.X() <<","<<uvV.Y() << " ";
1442 if ( hintLenOK || avgLenOK )
1444 // create a point between the closest node and V
1447 double r = min( 0.5, ( hintLenOK ? segLenHint/nodeDist3D : avgSegLen2d/nodeDist2D ));
1448 // direction from V to closet node in 2D
1449 gp_Dir2d v2n( helper.ApplyIn2D(surf, uvP, uvV, gp_XY_Subtracted, false ));
1451 uvP = vData.uv + r * nodeDist2D * v2n.XY();
1452 gp_Pnt P = surf->Value( uvP.X(), uvP.Y() ).Transformed( loc );
1454 netgen::MeshPoint mp( netgen::Point<3> (P.X(), P.Y(), P.Z()));
1455 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1456 segEnd2 = ngMesh.GetNP();
1457 //cout << "Middle " << r << " uv " << uvP.X() << "," << uvP.Y() << "( " << ngMesh.Point(segEnd2).X()<<","<<ngMesh.Point(segEnd2).Y()<<","<<ngMesh.Point(segEnd2).Z()<<" )"<< endl;
1458 SMDS_MeshNode * nP = helper.AddNode(P.X(), P.Y(), P.Z());
1459 nodeVec.push_back( nP );
1461 //else cout << "at Node " << " uv " << uvP.X() << "," << uvP.Y() << endl;
1464 netgen::Segment seg;
1466 if ( segEnd1 > segEnd2 ) swap( segEnd1, segEnd2 ), swap( uvV, uvP );
1467 seg[0] = segEnd1; // ng node id
1468 seg[1] = segEnd2; // ng node id
1469 seg.edgenr = ngMesh.GetNSeg() + 1;// segment id
1472 seg.epgeominfo[ 0 ].dist = 0; // param on curve
1473 seg.epgeominfo[ 0 ].u = uvV.X();
1474 seg.epgeominfo[ 0 ].v = uvV.Y();
1475 seg.epgeominfo[ 1 ].dist = 1; // param on curve
1476 seg.epgeominfo[ 1 ].u = uvP.X();
1477 seg.epgeominfo[ 1 ].v = uvP.Y();
1479 // seg.epgeominfo[ 0 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1480 // seg.epgeominfo[ 1 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1482 ngMesh.AddSegment (seg);
1484 // add reverse segment
1485 swap (seg[0], seg[1]);
1486 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1487 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1488 ngMesh.AddSegment (seg);
1494 //================================================================================
1496 * \brief Make netgen take internal vertices in solids into account by adding
1497 * faces including internal vertices
1499 * This function works in supposition that 2D mesh is already computed in ngMesh
1501 //================================================================================
1503 void NETGENPlugin_Mesher::AddIntVerticesInSolids(const netgen::OCCGeometry& occgeom,
1504 netgen::Mesh& ngMesh,
1505 vector<const SMDS_MeshNode*>& nodeVec,
1506 NETGENPlugin_Internals& internalShapes)
1508 #ifdef DUMP_TRIANGLES_SCRIPT
1509 // create a python script making a mesh containing triangles added for internal vertices
1510 ofstream py(DUMP_TRIANGLES_SCRIPT);
1511 py << "import SMESH"<< endl
1512 << "from salome.smesh import smeshBuilder"<<endl
1513 << "smesh = smeshBuilder.New(salome.myStudy)"<<endl
1514 << "m = smesh.Mesh(name='triangles')" << endl;
1516 if ((int) nodeVec.size() < ngMesh.GetNP() )
1517 nodeVec.resize( ngMesh.GetNP(), 0 );
1519 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1520 SMESH_MesherHelper helper( internalShapes.getMesh() );
1522 const map<int,list<int> >& so2Vert = internalShapes.getSolidsWithVertices();
1523 map<int,list<int> >::const_iterator s2v = so2Vert.begin();
1524 for ( ; s2v != so2Vert.end(); ++s2v )
1526 const TopoDS_Shape& solid = meshDS->IndexToShape( s2v->first );
1527 if ( solid.IsNull() ) continue;
1528 int solidNgID = occgeom.somap.FindIndex (solid);
1529 if ( solidNgID < 0 && !occgeom.somap.IsEmpty() ) continue;
1531 helper.SetSubShape( solid );
1532 helper.SetElementsOnShape( true );
1534 // find ng indices of faces within the solid
1536 for (TopExp_Explorer fExp(solid, TopAbs_FACE); fExp.More(); fExp.Next() )
1537 ngFaceIds.insert( occgeom.fmap.FindIndex( fExp.Current() ));
1538 if ( ngFaceIds.size() == 1 && *ngFaceIds.begin() == 0 )
1539 ngFaceIds.insert( 1 );
1541 // Get data of internal vertices and add them to ngMesh
1543 multimap< double, TIntVSoData > dist2VData; // sort vertices by distance from ng faces
1545 int i, nbFaceInit = ngMesh.GetNSE();
1547 // boundary characteristics
1548 double totSegLen = 0;
1551 const list<int>& iVertices = s2v->second;
1552 list<int>::const_iterator iv = iVertices.begin();
1553 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1556 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1558 // get node on vertex
1559 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1562 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1563 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1564 nV = SMESH_Algo::VertexNode( V, meshDS );
1565 if ( !nV ) continue;
1568 netgen::MeshPoint mpV( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1569 ngMesh.AddPoint ( mpV, 1, netgen::FIXEDPOINT );
1570 vData.ngId = ngMesh.GetNP();
1571 nodeVec.push_back( nV );
1573 // loop on all 2d elements to find the one closest to vertex and to count
1574 // average segment length
1575 double closeDist2 = numeric_limits<double>::max(), avgDist2;
1576 for (i = 1; i <= ngMesh.GetNSE(); ++i)
1578 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1579 if ( !ngFaceIds.count( elem.GetIndex() )) continue;
1581 multimap< double, int> dist2nID; // sort nodes of element by distance from V
1582 for ( int j = 0; j < elem.GetNP(); ++j)
1584 netgen::MeshPoint mp = ngMesh.Point( elem[j] );
1585 double d2 = dist2( mpV, mp );
1586 dist2nID.insert( make_pair( d2, elem[j] ));
1587 avgDist2 += d2 / elem.GetNP();
1589 totNbSeg++, totSegLen+= sqrt( dist2( mp, ngMesh.Point( elem[(j+1)%elem.GetNP()])));
1591 double dist = dist2nID.begin()->first; //avgDist2;
1592 if ( dist < closeDist2 )
1593 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= dist;
1595 dist2VData.insert( make_pair( closeDist2, vData ));
1598 if ( totNbSeg == 0 ) break;
1599 double avgSegLen = totSegLen / totNbSeg;
1601 // Loop on vertices to add triangles
1603 multimap< double, TIntVSoData >::iterator dist_vData = dist2VData.begin();
1604 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1606 double closeDist2 = dist_vData->first;
1607 TIntVSoData & vData = dist_vData->second;
1609 const netgen::MeshPoint& mpV = ngMesh.Point( vData.ngId );
1611 // try to find more close face among ones added for internal vertices
1612 for (i = nbFaceInit+1; i <= ngMesh.GetNSE(); ++i)
1614 double avgDist2 = 0;
1615 multimap< double, int> dist2nID;
1616 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1617 for ( int j = 0; j < elem.GetNP(); ++j)
1619 double d = dist2( mpV, ngMesh.Point( elem[j] ));
1620 dist2nID.insert( make_pair( d, elem[j] ));
1621 avgDist2 += d / elem.GetNP();
1622 if ( avgDist2 < closeDist2 )
1623 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= avgDist2;
1626 // sort nodes of the closest face by angle with vector from V to the closest node
1627 const double tol = numeric_limits<double>::min();
1628 map< double, int > angle2ID;
1629 const netgen::Element2d& closeFace = ngMesh.SurfaceElement( vData.ngIdClose );
1630 netgen::MeshPoint mp[2];
1631 mp[0] = ngMesh.Point( vData.ngIdCloseN );
1632 gp_XYZ p1( NGPOINT_COORDS( mp[0] ));
1633 gp_XYZ pV( NGPOINT_COORDS( mpV ));
1634 gp_Vec v2p1( pV, p1 );
1635 double distN1 = v2p1.Magnitude();
1636 if ( distN1 <= tol ) continue;
1638 for ( int j = 0; j < closeFace.GetNP(); ++j)
1640 mp[1] = ngMesh.Point( closeFace[j] );
1641 gp_Vec v2p( pV, gp_Pnt( NGPOINT_COORDS( mp[1] )) );
1642 angle2ID.insert( make_pair( v2p1.Angle( v2p ), closeFace[j]));
1644 // get node with angle of 60 degrees or greater
1645 map< double, int >::iterator angle_id = angle2ID.lower_bound( 60. * M_PI / 180. );
1646 if ( angle_id == angle2ID.end() ) angle_id = --angle2ID.end();
1647 const double minAngle = 30. * M_PI / 180.;
1648 const double angle = angle_id->first;
1649 bool angleOK = ( angle > minAngle );
1651 // find points to create a triangle
1652 netgen::Element2d tri(3);
1654 tri[0] = vData.ngId;
1655 tri[1] = vData.ngIdCloseN; // to use the closest nodes
1656 tri[2] = angle_id->second; // to use the node with best angle
1658 // decide whether to use the closest node and the node with best angle or to create new ones
1659 for ( int isBestAngleN = 0; isBestAngleN < 2; ++isBestAngleN )
1661 bool createNew = !angleOK; //, distOK = true;
1663 int triInd = isBestAngleN ? 2 : 1;
1664 mp[isBestAngleN] = ngMesh.Point( tri[triInd] );
1669 double distN2 = sqrt( dist2( mpV, mp[isBestAngleN]));
1670 createNew = ( fabs( distN2 - distN1 ) > 0.25 * distN1 );
1672 else if ( angle < tol )
1674 v2p1.SetX( v2p1.X() + 1e-3 );
1680 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1681 bool avgLenOK = ( avgSegLen < 0.75 * distN1 );
1682 bool hintLenOK = ( segLenHint < 0.75 * distN1 );
1683 createNew = (createNew || avgLenOK || hintLenOK );
1684 // we create a new node not closer than 0.5 to the closest face
1685 // in order not to clash with other close face
1686 double r = min( 0.5, ( hintLenOK ? segLenHint : avgSegLen ) / distN1 );
1687 distFromV = r * distN1;
1691 // create a new point, between the node and the vertex if angleOK
1692 gp_XYZ p( NGPOINT_COORDS( mp[isBestAngleN] ));
1693 gp_Vec v2p( pV, p ); v2p.Normalize();
1694 if ( isBestAngleN && !angleOK )
1695 p = p1 + gp_Dir( v2p.XYZ() - v2p1.XYZ()).XYZ() * distN1 * 0.95;
1697 p = pV + v2p.XYZ() * distFromV;
1699 if ( !isBestAngleN ) p1 = p, distN1 = distFromV;
1701 mp[isBestAngleN].SetPoint( netgen::Point<3> (p.X(), p.Y(), p.Z()));
1702 ngMesh.AddPoint ( mp[isBestAngleN], 1, netgen::SURFACEPOINT );
1703 tri[triInd] = ngMesh.GetNP();
1704 nodeVec.push_back( helper.AddNode( p.X(), p.Y(), p.Z()) );
1707 ngMesh.AddSurfaceElement (tri);
1708 swap( tri[1], tri[2] );
1709 ngMesh.AddSurfaceElement (tri);
1711 #ifdef DUMP_TRIANGLES_SCRIPT
1712 py << "n1 = m.AddNode( "<< mpV(0)<<", "<< mpV(1)<<", "<< mpV(2)<<") "<< endl
1713 << "n2 = m.AddNode( "<< mp[0](0)<<", "<< mp[0](1)<<", "<< mp[0](2)<<") "<< endl
1714 << "n3 = m.AddNode( "<< mp[1](0)<<", "<< mp[1](1)<<", "<< mp[1](2)<<" )" << endl
1715 << "m.AddFace([n1,n2,n3])" << endl;
1717 } // loop on internal vertices of a solid
1719 } // loop on solids with internal vertices
1722 //================================================================================
1724 * \brief Fill netgen mesh with segments of a FACE
1725 * \param ngMesh - netgen mesh
1726 * \param geom - container of OCCT geometry to mesh
1727 * \param wires - data of nodes on FACE boundary
1728 * \param helper - mesher helper holding the FACE
1729 * \param nodeVec - vector of nodes in which node index == netgen ID
1730 * \retval SMESH_ComputeErrorPtr - error description
1732 //================================================================================
1734 SMESH_ComputeErrorPtr
1735 NETGENPlugin_Mesher::AddSegmentsToMesh(netgen::Mesh& ngMesh,
1736 netgen::OCCGeometry& geom,
1737 const TSideVector& wires,
1738 SMESH_MesherHelper& helper,
1739 vector< const SMDS_MeshNode* > & nodeVec,
1740 const bool overrideMinH)
1742 // ----------------------------
1743 // Check wires and count nodes
1744 // ----------------------------
1746 for ( size_t iW = 0; iW < wires.size(); ++iW )
1748 StdMeshers_FaceSidePtr wire = wires[ iW ];
1749 if ( wire->MissVertexNode() )
1751 // Commented for issue 0020960. It worked for the case, let's wait for case where it doesn't.
1752 // It seems that there is no reason for this limitation
1754 // (new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH, "Missing nodes on vertices"));
1756 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1757 if ((int) uvPtVec.size() != wire->NbPoints() )
1758 return SMESH_ComputeError::New(COMPERR_BAD_INPUT_MESH,
1759 SMESH_Comment("Unexpected nb of points on wire ") << iW
1760 << ": " << uvPtVec.size()<<" != "<<wire->NbPoints());
1761 nbNodes += wire->NbPoints();
1763 nodeVec.reserve( nodeVec.size() + nbNodes + 1 );
1764 if ( nodeVec.empty() )
1765 nodeVec.push_back( 0 );
1767 // -----------------
1769 // -----------------
1771 const bool wasNgMeshEmpty = ( ngMesh.GetNP() < 1 ); /* true => this method is called by
1772 NETGENPlugin_NETGEN_2D_ONLY */
1774 // map for nodes on vertices since they can be shared between wires
1775 // ( issue 0020676, face_int_box.brep) and nodes built by NETGEN
1776 map<const SMDS_MeshNode*, int > node2ngID;
1777 if ( !wasNgMeshEmpty ) // fill node2ngID with nodes built by NETGEN
1779 set< int > subIDs; // ids of sub-shapes of the FACE
1780 for ( size_t iW = 0; iW < wires.size(); ++iW )
1782 StdMeshers_FaceSidePtr wire = wires[ iW ];
1783 for ( int iE = 0, nbE = wire->NbEdges(); iE < nbE; ++iE )
1785 subIDs.insert( wire->EdgeID( iE ));
1786 subIDs.insert( helper.GetMeshDS()->ShapeToIndex( wire->FirstVertex( iE )));
1789 for ( size_t ngID = 1; ngID < nodeVec.size(); ++ngID )
1790 if ( subIDs.count( nodeVec[ngID]->getshapeId() ))
1791 node2ngID.insert( make_pair( nodeVec[ngID], ngID ));
1794 const int solidID = 0, faceID = geom.fmap.FindIndex( helper.GetSubShape() );
1795 if ( ngMesh.GetNFD() < 1 )
1796 ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(faceID, solidID, solidID, 0));
1798 for ( size_t iW = 0; iW < wires.size(); ++iW )
1800 StdMeshers_FaceSidePtr wire = wires[ iW ];
1801 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1802 const int nbSegments = wire->NbPoints() - 1;
1804 // assure the 1st node to be in node2ngID, which is needed to correctly
1805 // "close chain of segments" (see below) in case if the 1st node is not
1806 // onVertex because it is on a Viscous layer
1807 node2ngID.insert( make_pair( uvPtVec[ 0 ].node, ngMesh.GetNP() + 1 ));
1809 // compute length of every segment
1810 vector<double> segLen( nbSegments );
1811 for ( int i = 0; i < nbSegments; ++i )
1812 segLen[i] = SMESH_TNodeXYZ( uvPtVec[ i ].node ).Distance( uvPtVec[ i+1 ].node );
1814 int edgeID = 1, posID = -2;
1815 bool isInternalWire = false;
1816 double vertexNormPar = 0;
1817 //const int prevNbNGSeg = ngMesh.GetNSeg();
1818 for ( int i = 0; i < nbSegments; ++i ) // loop on segments
1820 // Add the first point of a segment
1822 const SMDS_MeshNode * n = uvPtVec[ i ].node;
1823 const int posShapeID = n->getshapeId();
1824 bool onVertex = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX );
1825 bool onEdge = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_EDGE );
1827 // skip nodes on degenerated edges
1828 if ( helper.IsDegenShape( posShapeID ) &&
1829 helper.IsDegenShape( uvPtVec[ i+1 ].node->getshapeId() ))
1832 int ngID1 = ngMesh.GetNP() + 1, ngID2 = ngID1+1;
1833 if ( onVertex || ( !wasNgMeshEmpty && onEdge ) || helper.IsRealSeam( posShapeID ))
1834 ngID1 = node2ngID.insert( make_pair( n, ngID1 )).first->second;
1835 if ( ngID1 > ngMesh.GetNP() )
1837 netgen::MeshPoint mp( netgen::Point<3> (n->X(), n->Y(), n->Z()) );
1838 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1839 nodeVec.push_back( n );
1841 else // n is in ngMesh already, and ngID2 in prev segment is wrong
1843 ngID2 = ngMesh.GetNP() + 1;
1844 if ( i > 0 ) // prev segment belongs to same wire
1846 netgen::Segment& prevSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
1853 netgen::Segment seg;
1855 seg[0] = ngID1; // ng node id
1856 seg[1] = ngID2; // ng node id
1857 seg.edgenr = ngMesh.GetNSeg() + 1; // ng segment id
1858 seg.si = faceID; // = geom.fmap.FindIndex (face);
1860 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1862 const UVPtStruct& pnt = uvPtVec[ i + iEnd ];
1864 seg.epgeominfo[ iEnd ].dist = pnt.param; // param on curve
1865 seg.epgeominfo[ iEnd ].u = pnt.u;
1866 seg.epgeominfo[ iEnd ].v = pnt.v;
1868 // find out edge id and node parameter on edge
1869 onVertex = ( pnt.normParam + 1e-10 > vertexNormPar );
1870 if ( onVertex || posShapeID != posID )
1873 double normParam = pnt.normParam;
1875 normParam = 0.5 * ( uvPtVec[ i ].normParam + uvPtVec[ i+1 ].normParam );
1876 int edgeIndexInWire = wire->EdgeIndex( normParam );
1877 vertexNormPar = wire->LastParameter( edgeIndexInWire );
1878 const TopoDS_Edge& edge = wire->Edge( edgeIndexInWire );
1879 edgeID = geom.emap.FindIndex( edge );
1881 isInternalWire = ( edge.Orientation() == TopAbs_INTERNAL );
1882 // if ( onVertex ) // param on curve is different on each of two edges
1883 // seg.epgeominfo[ iEnd ].dist = helper.GetNodeU( edge, pnt.node );
1885 seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
1888 ngMesh.AddSegment (seg);
1890 // restrict size of elements near the segment
1891 SMESH_TNodeXYZ np1( n ), np2( uvPtVec[ i+1 ].node );
1892 // get an average size of adjacent segments to avoid sharp change of
1893 // element size (regression on issue 0020452, note 0010898)
1894 int iPrev = SMESH_MesherHelper::WrapIndex( i-1, nbSegments );
1895 int iNext = SMESH_MesherHelper::WrapIndex( i+1, nbSegments );
1896 double sumH = segLen[ iPrev ] + segLen[ i ] + segLen[ iNext ];
1897 int nbSeg = ( int( segLen[ iPrev ] > sumH / 100.) +
1898 int( segLen[ i ] > sumH / 100.) +
1899 int( segLen[ iNext ] > sumH / 100.));
1901 RestrictLocalSize( ngMesh, 0.5*(np1+np2), sumH / nbSeg, overrideMinH );
1903 if ( isInternalWire )
1905 swap (seg[0], seg[1]);
1906 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1907 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1908 ngMesh.AddSegment (seg);
1910 } // loop on segments on a wire
1912 // close chain of segments
1913 if ( nbSegments > 0 )
1915 netgen::Segment& lastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() - int( isInternalWire));
1916 const SMDS_MeshNode * lastNode = uvPtVec.back().node;
1917 lastSeg[1] = node2ngID.insert( make_pair( lastNode, lastSeg[1] )).first->second;
1918 if ( lastSeg[1] > ngMesh.GetNP() )
1920 netgen::MeshPoint mp( netgen::Point<3> (lastNode->X(), lastNode->Y(), lastNode->Z()) );
1921 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1922 nodeVec.push_back( lastNode );
1924 if ( isInternalWire )
1926 netgen::Segment& realLastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
1927 realLastSeg[0] = lastSeg[1];
1931 #ifdef DUMP_SEGMENTS
1932 cout << "BEGIN WIRE " << iW << endl;
1933 for ( int i = prevNbNGSeg+1; i <= ngMesh.GetNSeg(); ++i )
1935 netgen::Segment& seg = ngMesh.LineSegment( i );
1937 netgen::Segment& prevSeg = ngMesh.LineSegment( i-1 );
1938 if ( seg[0] == prevSeg[1] && seg[1] == prevSeg[0] )
1940 cout << "Segment: " << seg.edgenr << endl << "\tis REVRESE of the previous one" << endl;
1944 cout << "Segment: " << seg.edgenr << endl
1945 << "\tp1: " << seg[0] << " n" << nodeVec[ seg[0]]->GetID() << endl
1946 << "\tp2: " << seg[1] << " n" << nodeVec[ seg[1]]->GetID() << endl
1947 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
1948 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
1949 << "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
1950 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
1951 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl
1952 << "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
1954 cout << "--END WIRE " << iW << endl;
1957 } // loop on WIREs of a FACE
1959 // add a segment instead of an internal vertex
1960 if ( wasNgMeshEmpty )
1962 NETGENPlugin_Internals intShapes( *helper.GetMesh(), helper.GetSubShape(), /*is3D=*/false );
1963 AddIntVerticesInFaces( geom, ngMesh, nodeVec, intShapes );
1965 ngMesh.CalcSurfacesOfNode();
1970 //================================================================================
1972 * \brief Fill SMESH mesh according to contents of netgen mesh
1973 * \param occgeo - container of OCCT geometry to mesh
1974 * \param ngMesh - netgen mesh
1975 * \param initState - bn of entities in netgen mesh before computing
1976 * \param sMesh - SMESH mesh to fill in
1977 * \param nodeVec - vector of nodes in which node index == netgen ID
1978 * \param comment - returns problem description
1979 * \param quadHelper - holder of medium nodes of sub-meshes
1980 * \retval int - error
1982 //================================================================================
1984 int NETGENPlugin_Mesher::FillSMesh(const netgen::OCCGeometry& occgeo,
1985 netgen::Mesh& ngMesh,
1986 const NETGENPlugin_ngMeshInfo& initState,
1988 std::vector<const SMDS_MeshNode*>& nodeVec,
1989 SMESH_Comment& comment,
1990 SMESH_MesherHelper* quadHelper)
1992 int nbNod = ngMesh.GetNP();
1993 int nbSeg = ngMesh.GetNSeg();
1994 int nbFac = ngMesh.GetNSE();
1995 int nbVol = ngMesh.GetNE();
1997 SMESHDS_Mesh* meshDS = sMesh.GetMeshDS();
1999 // quadHelper is used for either
2000 // 1) making quadratic elements when a lower dimention mesh is loaded
2001 // to SMESH before convertion to quadratic by NETGEN
2002 // 2) sewing of quadratic elements with quadratic elements of sub-meshes
2003 if ( quadHelper && !quadHelper->GetIsQuadratic() && quadHelper->GetTLinkNodeMap().empty() )
2006 // -------------------------------------
2007 // Create and insert nodes into nodeVec
2008 // -------------------------------------
2010 nodeVec.resize( nbNod + 1 );
2011 int i, nbInitNod = initState._nbNodes;
2012 for (i = nbInitNod+1; i <= nbNod; ++i )
2014 const netgen::MeshPoint& ngPoint = ngMesh.Point(i);
2015 SMDS_MeshNode* node = NULL;
2016 TopoDS_Vertex aVert;
2017 // First, netgen creates nodes on vertices in occgeo.vmap,
2018 // so node index corresponds to vertex index
2019 // but (issue 0020776) netgen does not create nodes with equal coordinates
2020 if ( i-nbInitNod <= occgeo.vmap.Extent() )
2022 gp_Pnt p ( NGPOINT_COORDS(ngPoint) );
2023 for (int iV = i-nbInitNod; aVert.IsNull() && iV <= occgeo.vmap.Extent(); ++iV)
2025 aVert = TopoDS::Vertex( occgeo.vmap( iV ));
2026 gp_Pnt pV = BRep_Tool::Pnt( aVert );
2027 if ( p.SquareDistance( pV ) > 1e-20 )
2030 node = const_cast<SMDS_MeshNode*>( SMESH_Algo::VertexNode( aVert, meshDS ));
2033 if (!node) // node not found on vertex
2035 node = meshDS->AddNode( NGPOINT_COORDS( ngPoint ));
2036 if (!aVert.IsNull())
2037 meshDS->SetNodeOnVertex(node, aVert);
2042 // -------------------------------------------
2043 // Create mesh segments along geometric edges
2044 // -------------------------------------------
2046 int nbInitSeg = initState._nbSegments;
2047 for (i = nbInitSeg+1; i <= nbSeg; ++i )
2049 const netgen::Segment& seg = ngMesh.LineSegment(i);
2051 int pinds[3] = { seg.pnums[0], seg.pnums[1], seg.pnums[2] };
2054 for (int j=0; j < 3; ++j)
2056 int pind = pinds[j];
2057 if (pind <= 0 || !nodeVec_ACCESS(pind))
2065 int aGeomEdgeInd = seg.epgeominfo[j].edgenr;
2066 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
2067 aEdge = TopoDS::Edge(occgeo.emap(aGeomEdgeInd));
2069 param = seg.epgeominfo[j].dist;
2072 else // middle point
2074 param = param2 * 0.5;
2076 if (!aEdge.IsNull() && nodeVec_ACCESS(pind)->getshapeId() < 1)
2078 meshDS->SetNodeOnEdge(nodeVec_ACCESS(pind), aEdge, param);
2083 SMDS_MeshEdge* edge = 0;
2084 if (nbp == 2) // second order ?
2086 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1])))
2088 if ( quadHelper ) // final mesh must be quadratic
2089 edge = quadHelper->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
2091 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
2095 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
2096 nodeVec_ACCESS(pinds[2])))
2098 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
2099 nodeVec_ACCESS(pinds[2]));
2103 if ( comment.empty() ) comment << "Cannot create a mesh edge";
2104 MESSAGE("Cannot create a mesh edge");
2105 nbSeg = nbFac = nbVol = 0;
2108 if ( !aEdge.IsNull() && edge->getshapeId() < 1 )
2109 meshDS->SetMeshElementOnShape(edge, aEdge);
2111 else if ( comment.empty() )
2113 comment << "Invalid netgen segment #" << i;
2117 // ----------------------------------------
2118 // Create mesh faces along geometric faces
2119 // ----------------------------------------
2121 int nbInitFac = initState._nbFaces;
2122 int quadFaceID = ngMesh.GetNFD() + 1;
2123 if ( nbInitFac < nbFac )
2124 // add a faces descriptor to exclude qudrangle elements generated by NETGEN
2125 // from computation of 3D mesh
2126 ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(quadFaceID, /*solid1=*/0, /*solid2=*/0, 0));
2128 vector<const SMDS_MeshNode*> nodes;
2129 for (i = nbInitFac+1; i <= nbFac; ++i )
2131 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
2132 int aGeomFaceInd = elem.GetIndex();
2134 if (aGeomFaceInd > 0 && aGeomFaceInd <= occgeo.fmap.Extent())
2135 aFace = TopoDS::Face(occgeo.fmap(aGeomFaceInd));
2137 for (int j=1; j <= elem.GetNP(); ++j)
2139 int pind = elem.PNum(j);
2140 if ( pind < 1 || pind >= (int) nodeVec.size() )
2142 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind))
2144 nodes.push_back( node );
2145 if (!aFace.IsNull() && node->getshapeId() < 1)
2147 const netgen::PointGeomInfo& pgi = elem.GeomInfoPi(j);
2148 meshDS->SetNodeOnFace(node, aFace, pgi.u, pgi.v);
2152 if ((int) nodes.size() != elem.GetNP() )
2154 if ( comment.empty() )
2155 comment << "Invalid netgen 2d element #" << i;
2156 continue; // bad node ids
2158 SMDS_MeshFace* face = NULL;
2159 switch (elem.GetType())
2162 if ( quadHelper ) // final mesh must be quadratic
2163 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2]);
2165 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2]);
2168 if ( quadHelper ) // final mesh must be quadratic
2169 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2171 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2172 // exclude qudrangle elements from computation of 3D mesh
2173 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2176 nodes[5] = mediumNode( nodes[0],nodes[1],nodes[5], quadHelper );
2177 nodes[3] = mediumNode( nodes[1],nodes[2],nodes[3], quadHelper );
2178 nodes[4] = mediumNode( nodes[2],nodes[0],nodes[4], quadHelper );
2179 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[5],nodes[3],nodes[4]);
2182 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2183 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2184 nodes[5] = mediumNode( nodes[2],nodes[3],nodes[5], quadHelper );
2185 nodes[6] = mediumNode( nodes[3],nodes[0],nodes[6], quadHelper );
2186 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3],
2187 nodes[4],nodes[7],nodes[5],nodes[6]);
2188 // exclude qudrangle elements from computation of 3D mesh
2189 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2192 MESSAGE("NETGEN created a face of unexpected type, ignoring");
2197 if ( comment.empty() ) comment << "Cannot create a mesh face";
2198 MESSAGE("Cannot create a mesh face");
2199 nbSeg = nbFac = nbVol = 0;
2202 if (!aFace.IsNull())
2203 meshDS->SetMeshElementOnShape(face, aFace);
2206 // ------------------
2207 // Create tetrahedra
2208 // ------------------
2210 for (i = 1; i <= nbVol; ++i)
2212 const netgen::Element& elem = ngMesh.VolumeElement(i);
2213 int aSolidInd = elem.GetIndex();
2214 TopoDS_Solid aSolid;
2215 if (aSolidInd > 0 && aSolidInd <= occgeo.somap.Extent())
2216 aSolid = TopoDS::Solid(occgeo.somap(aSolidInd));
2218 for (int j=1; j <= elem.GetNP(); ++j)
2220 int pind = elem.PNum(j);
2221 if ( pind < 1 || pind >= (int)nodeVec.size() )
2223 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind) )
2225 nodes.push_back(node);
2226 if ( !aSolid.IsNull() && node->getshapeId() < 1 )
2227 meshDS->SetNodeInVolume(node, aSolid);
2230 if ((int) nodes.size() != elem.GetNP() )
2232 if ( comment.empty() )
2233 comment << "Invalid netgen 3d element #" << i;
2236 SMDS_MeshVolume* vol = NULL;
2237 switch (elem.GetType())
2240 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3]);
2243 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2244 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2245 nodes[5] = mediumNode( nodes[2],nodes[0],nodes[5], quadHelper );
2246 nodes[6] = mediumNode( nodes[0],nodes[3],nodes[6], quadHelper );
2247 nodes[8] = mediumNode( nodes[1],nodes[3],nodes[8], quadHelper );
2248 nodes[9] = mediumNode( nodes[2],nodes[3],nodes[9], quadHelper );
2249 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3],
2250 nodes[4],nodes[7],nodes[5],nodes[6],nodes[8],nodes[9]);
2253 MESSAGE("NETGEN created a volume of unexpected type, ignoring");
2258 if ( comment.empty() ) comment << "Cannot create a mesh volume";
2259 MESSAGE("Cannot create a mesh volume");
2260 nbSeg = nbFac = nbVol = 0;
2263 if (!aSolid.IsNull())
2264 meshDS->SetMeshElementOnShape(vol, aSolid);
2266 return comment.empty() ? 0 : 1;
2271 //================================================================================
2273 * \brief Restrict size of elements on the given edge
2275 //================================================================================
2277 void setLocalSize(const TopoDS_Edge& edge,
2281 if ( size <= std::numeric_limits<double>::min() )
2283 Standard_Real u1, u2;
2284 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, u1, u2);
2285 if ( curve.IsNull() )
2287 TopoDS_Iterator vIt( edge );
2288 if ( !vIt.More() ) return;
2289 gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex( vIt.Value() ));
2290 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
2294 const int nb = (int)( 1.5 * SMESH_Algo::EdgeLength( edge ) / size );
2295 Standard_Real delta = (u2-u1)/nb;
2296 for(int i=0; i<nb; i++)
2298 Standard_Real u = u1 + delta*i;
2299 gp_Pnt p = curve->Value(u);
2300 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
2301 netgen::Point3d pi(p.X(), p.Y(), p.Z());
2302 double resultSize = mesh.GetH(pi);
2303 if ( resultSize - size > 0.1*size )
2304 // netgen does restriction iff oldH/newH > 1.2 (localh.cpp:136)
2305 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), resultSize/1.201 );
2310 //================================================================================
2312 * \brief Convert error into text
2314 //================================================================================
2316 std::string text(int err)
2321 SMESH_Comment("Error in netgen::OCCGenerateMesh() at ") << netgen::multithread.task;
2324 //================================================================================
2326 * \brief Convert exception into text
2328 //================================================================================
2330 std::string text(Standard_Failure& ex)
2332 SMESH_Comment str("Exception in netgen::OCCGenerateMesh()");
2333 str << " at " << netgen::multithread.task
2334 << ": " << ex.DynamicType()->Name();
2335 if ( ex.GetMessageString() && strlen( ex.GetMessageString() ))
2336 str << ": " << ex.GetMessageString();
2339 //================================================================================
2341 * \brief Convert exception into text
2343 //================================================================================
2345 std::string text(netgen::NgException& ex)
2347 SMESH_Comment str("NgException");
2348 if ( strlen( netgen::multithread.task ) > 0 )
2349 str << " at " << netgen::multithread.task;
2350 str << ": " << ex.What();
2354 //================================================================================
2356 * \brief Looks for triangles lying on a SOLID
2358 //================================================================================
2360 bool hasBadElemOnSolid( const list<const SMDS_MeshElement*>& elems,
2361 SMESH_subMesh* solidSM )
2363 TopTools_IndexedMapOfShape solidSubs;
2364 TopExp::MapShapes( solidSM->GetSubShape(), solidSubs );
2365 SMESHDS_Mesh* mesh = solidSM->GetFather()->GetMeshDS();
2367 list<const SMDS_MeshElement*>::const_iterator e = elems.begin();
2368 for ( ; e != elems.end(); ++e )
2370 const SMDS_MeshElement* elem = *e;
2371 // if ( elem->GetType() != SMDSAbs_Face ) -- 23047
2373 int nbNodesOnSolid = 0, nbNodes = elem->NbNodes();
2374 SMDS_NodeIteratorPtr nIt = elem->nodeIterator();
2375 while ( nIt->more() )
2377 const SMDS_MeshNode* n = nIt->next();
2378 const TopoDS_Shape& s = mesh->IndexToShape( n->getshapeId() );
2379 nbNodesOnSolid += ( !s.IsNull() && solidSubs.Contains( s ));
2380 if ( nbNodesOnSolid > 2 ||
2381 nbNodesOnSolid == nbNodes)
2388 const double edgeMeshingTime = 0.001;
2389 const double faceMeshingTime = 0.019;
2390 const double edgeFaceMeshingTime = edgeMeshingTime + faceMeshingTime;
2391 const double faceOptimizTime = 0.06;
2392 const double voluMeshingTime = 0.15;
2393 const double volOptimizeTime = 0.77;
2396 //=============================================================================
2398 * Here we are going to use the NETGEN mesher
2400 //=============================================================================
2402 bool NETGENPlugin_Mesher::Compute()
2404 NETGENPlugin_NetgenLibWrapper ngLib;
2406 netgen::MeshingParameters& mparams = netgen::mparam;
2407 MESSAGE("Compute with:\n"
2408 " max size = " << mparams.maxh << "\n"
2409 " segments per edge = " << mparams.segmentsperedge);
2411 " growth rate = " << mparams.grading << "\n"
2412 " elements per radius = " << mparams.curvaturesafety << "\n"
2413 " second order = " << mparams.secondorder << "\n"
2414 " quad allowed = " << mparams.quad << "\n"
2415 " surface curvature = " << mparams.uselocalh << "\n"
2416 " fuse edges = " << netgen::merge_solids);
2418 SMESH_ComputeErrorPtr error = SMESH_ComputeError::New();
2419 SMESH_MesherHelper quadHelper( *_mesh );
2420 quadHelper.SetIsQuadratic( mparams.secondorder );
2422 static string debugFile = "/tmp/ngMesh.py"; /* to call toPython( _ngMesh, debugFile )
2423 while debugging netgen */
2424 // -------------------------
2425 // Prepare OCC geometry
2426 // -------------------------
2428 netgen::OCCGeometry occgeo;
2429 list< SMESH_subMesh* > meshedSM[3]; // for 0-2 dimensions
2430 NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
2431 PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
2434 _totalTime = edgeFaceMeshingTime;
2436 _totalTime += faceOptimizTime;
2438 _totalTime += voluMeshingTime + ( _optimize ? volOptimizeTime : 0 );
2439 double doneTime = 0;
2442 _curShapeIndex = -1;
2444 // -------------------------
2445 // Generate the mesh
2446 // -------------------------
2449 NETGENPlugin_ngMeshInfo initState; // it remembers size of ng mesh equal to size of Smesh
2451 SMESH_Comment comment;
2454 // vector of nodes in which node index == netgen ID
2455 vector< const SMDS_MeshNode* > nodeVec;
2463 // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
2464 mparams.uselocalh = false;
2465 mparams.grading = 0.8; // not limitited size growth
2467 if ( _simpleHyp->GetNumberOfSegments() )
2469 mparams.maxh = occgeo.boundingbox.Diam();
2472 mparams.maxh = _simpleHyp->GetLocalLength();
2475 if ( mparams.maxh == 0.0 )
2476 mparams.maxh = occgeo.boundingbox.Diam();
2477 if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
2478 mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
2480 // Local size on faces
2481 occgeo.face_maxh = mparams.maxh;
2483 // Let netgen create _ngMesh and calculate element size on not meshed shapes
2487 int startWith = netgen::MESHCONST_ANALYSE;
2488 int endWith = netgen::MESHCONST_ANALYSE;
2493 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2495 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2497 if(netgen::multithread.terminate)
2500 comment << text(err);
2502 catch (Standard_Failure& ex)
2504 comment << text(ex);
2506 err = 0; //- MESHCONST_ANALYSE isn't so important step
2509 ngLib.setMesh(( Ng_Mesh*) _ngMesh );
2511 _ngMesh->ClearFaceDescriptors(); // we make descriptors our-self
2515 // Pass 1D simple parameters to NETGEN
2516 // --------------------------------
2517 int nbSeg = _simpleHyp->GetNumberOfSegments();
2518 double segSize = _simpleHyp->GetLocalLength();
2519 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
2521 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
2523 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
2524 setLocalSize( e, segSize, *_ngMesh );
2527 else // if ( ! _simpleHyp )
2529 // Local size on vertices and edges
2530 // --------------------------------
2531 for(std::map<int,double>::const_iterator it=EdgeId2LocalSize.begin(); it!=EdgeId2LocalSize.end(); it++)
2533 int key = (*it).first;
2534 double hi = (*it).second;
2535 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2536 const TopoDS_Edge& e = TopoDS::Edge(shape);
2537 setLocalSize( e, hi, *_ngMesh );
2539 for(std::map<int,double>::const_iterator it=VertexId2LocalSize.begin(); it!=VertexId2LocalSize.end(); it++)
2541 int key = (*it).first;
2542 double hi = (*it).second;
2543 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2544 const TopoDS_Vertex& v = TopoDS::Vertex(shape);
2545 gp_Pnt p = BRep_Tool::Pnt(v);
2546 NETGENPlugin_Mesher::RestrictLocalSize( *_ngMesh, p.XYZ(), hi );
2548 for(map<int,double>::const_iterator it=FaceId2LocalSize.begin();
2549 it!=FaceId2LocalSize.end(); it++)
2551 int key = (*it).first;
2552 double val = (*it).second;
2553 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2554 int faceNgID = occgeo.fmap.FindIndex(shape);
2555 occgeo.SetFaceMaxH(faceNgID, val);
2556 for ( TopExp_Explorer edgeExp( shape, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
2557 setLocalSize( TopoDS::Edge( edgeExp.Current() ), val, *_ngMesh );
2561 // Precompute internal edges (issue 0020676) in order to
2562 // add mesh on them correctly (twice) to netgen mesh
2563 if ( !err && internals.hasInternalEdges() )
2565 // load internal shapes into OCCGeometry
2566 netgen::OCCGeometry intOccgeo;
2567 internals.getInternalEdges( intOccgeo.fmap, intOccgeo.emap, intOccgeo.vmap, meshedSM );
2568 intOccgeo.boundingbox = occgeo.boundingbox;
2569 intOccgeo.shape = occgeo.shape;
2570 intOccgeo.face_maxh.SetSize(intOccgeo.fmap.Extent());
2571 intOccgeo.face_maxh = netgen::mparam.maxh;
2572 netgen::Mesh *tmpNgMesh = NULL;
2576 // compute local H on internal shapes in the main mesh
2577 //OCCSetLocalMeshSize(intOccgeo, *_ngMesh); it deletes _ngMesh->localH
2579 // let netgen create a temporary mesh
2581 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2583 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2585 if(netgen::multithread.terminate)
2588 // copy LocalH from the main to temporary mesh
2589 initState.transferLocalH( _ngMesh, tmpNgMesh );
2591 // compute mesh on internal edges
2592 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2594 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2596 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2598 comment << text(err);
2600 catch (Standard_Failure& ex)
2602 comment << text(ex);
2605 initState.restoreLocalH( tmpNgMesh );
2607 // fill SMESH by netgen mesh
2608 vector< const SMDS_MeshNode* > tmpNodeVec;
2609 FillSMesh( intOccgeo, *tmpNgMesh, initState, *_mesh, tmpNodeVec, comment );
2610 err = ( err || !comment.empty() );
2612 nglib::Ng_DeleteMesh((nglib::Ng_Mesh*)tmpNgMesh);
2615 // Fill _ngMesh with nodes and segments of computed submeshes
2618 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_0D ]) &&
2619 FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_1D ], &quadHelper));
2621 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2626 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2631 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2633 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2635 if(netgen::multithread.terminate)
2638 comment << text(err);
2640 catch (Standard_Failure& ex)
2642 comment << text(ex);
2647 _ticTime = ( doneTime += edgeMeshingTime ) / _totalTime / _progressTic;
2649 mparams.uselocalh = true; // restore as it is used at surface optimization
2651 // ---------------------
2652 // compute surface mesh
2653 // ---------------------
2656 // Pass 2D simple parameters to NETGEN
2658 if ( double area = _simpleHyp->GetMaxElementArea() ) {
2660 mparams.maxh = sqrt(2. * area/sqrt(3.0));
2661 mparams.grading = 0.4; // moderate size growth
2664 // length from edges
2665 if ( _ngMesh->GetNSeg() ) {
2666 double edgeLength = 0;
2667 TopTools_MapOfShape visitedEdges;
2668 for ( TopExp_Explorer exp( _shape, TopAbs_EDGE ); exp.More(); exp.Next() )
2669 if( visitedEdges.Add(exp.Current()) )
2670 edgeLength += SMESH_Algo::EdgeLength( TopoDS::Edge( exp.Current() ));
2671 // we have to multiply length by 2 since for each TopoDS_Edge there
2672 // are double set of NETGEN edges, in other words, we have to
2673 // divide _ngMesh->GetNSeg() by 2.
2674 mparams.maxh = 2*edgeLength / _ngMesh->GetNSeg();
2677 mparams.maxh = 1000;
2679 mparams.grading = 0.2; // slow size growth
2681 mparams.quad = _simpleHyp->GetAllowQuadrangles();
2682 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2683 _ngMesh->SetGlobalH (mparams.maxh);
2684 netgen::Box<3> bb = occgeo.GetBoundingBox();
2685 bb.Increase (bb.Diam()/20);
2686 _ngMesh->SetLocalH (bb.PMin(), bb.PMax(), mparams.grading);
2689 // Care of vertices internal in faces (issue 0020676)
2690 if ( internals.hasInternalVertexInFace() )
2692 // store computed segments in SMESH in order not to create SMESH
2693 // edges for ng segments added by AddIntVerticesInFaces()
2694 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2695 // add segments to faces with internal vertices
2696 AddIntVerticesInFaces( occgeo, *_ngMesh, nodeVec, internals );
2697 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2700 // Build viscous layers
2701 if ( _isViscousLayers2D )
2703 if ( !internals.hasInternalVertexInFace() ) {
2704 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2705 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2707 SMESH_ProxyMesh::Ptr viscousMesh;
2708 SMESH_MesherHelper helper( *_mesh );
2709 for ( int faceID = 1; faceID <= occgeo.fmap.Extent(); ++faceID )
2711 const TopoDS_Face& F = TopoDS::Face( occgeo.fmap( faceID ));
2712 viscousMesh = StdMeshers_ViscousLayers2D::Compute( *_mesh, F );
2715 // exclude from computation ng segments built on EDGEs of F
2716 for (int i = 1; i <= _ngMesh->GetNSeg(); i++)
2718 netgen::Segment & seg = _ngMesh->LineSegment(i);
2719 if (seg.si == faceID)
2722 // add new segments to _ngMesh instead of excluded ones
2723 helper.SetSubShape( F );
2725 StdMeshers_FaceSide::GetFaceWires( F, *_mesh, /*skipMediumNodes=*/true,
2726 error, viscousMesh );
2727 error = AddSegmentsToMesh( *_ngMesh, occgeo, wires, helper, nodeVec );
2729 if ( !error ) error = SMESH_ComputeError::New();
2731 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2734 // Let netgen compute 2D mesh
2735 startWith = netgen::MESHCONST_MESHSURFACE;
2736 endWith = _optimize ? netgen::MESHCONST_OPTSURFACE : netgen::MESHCONST_MESHSURFACE;
2741 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2743 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2745 if(netgen::multithread.terminate)
2748 comment << text (err);
2750 catch (Standard_Failure& ex)
2752 comment << text(ex);
2753 //err = 1; -- try to make volumes anyway
2755 catch (netgen::NgException exc)
2757 comment << text(exc);
2758 //err = 1; -- try to make volumes anyway
2763 doneTime += faceMeshingTime + ( _optimize ? faceOptimizTime : 0 );
2764 _ticTime = doneTime / _totalTime / _progressTic;
2766 // ---------------------
2767 // generate volume mesh
2768 // ---------------------
2769 // Fill _ngMesh with nodes and faces of computed 2D submeshes
2770 if ( !err && _isVolume && ( !meshedSM[ MeshDim_2D ].empty() || mparams.quad ))
2772 // load SMESH with computed segments and faces
2773 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2775 // compute pyramids on quadrangles
2776 SMESH_ProxyMesh::Ptr proxyMesh;
2777 if ( _mesh->NbQuadrangles() > 0 )
2778 for ( int iS = 1; iS <= occgeo.somap.Extent(); ++iS )
2780 StdMeshers_QuadToTriaAdaptor* Adaptor = new StdMeshers_QuadToTriaAdaptor;
2781 proxyMesh.reset( Adaptor );
2783 int nbPyrams = _mesh->NbPyramids();
2784 Adaptor->Compute( *_mesh, occgeo.somap(iS) );
2785 if ( nbPyrams != _mesh->NbPyramids() )
2787 list< SMESH_subMesh* > quadFaceSM;
2788 for (TopExp_Explorer face(occgeo.somap(iS), TopAbs_FACE); face.More(); face.Next())
2789 if ( Adaptor->GetProxySubMesh( face.Current() ))
2791 quadFaceSM.push_back( _mesh->GetSubMesh( face.Current() ));
2792 meshedSM[ MeshDim_2D ].remove( quadFaceSM.back() );
2794 FillNgMesh(occgeo, *_ngMesh, nodeVec, quadFaceSM, &quadHelper, proxyMesh);
2797 // fill _ngMesh with faces of sub-meshes
2798 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_2D ], &quadHelper));
2799 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2800 //toPython( _ngMesh, "/tmp/ngPython.py");
2802 if (!err && _isVolume)
2804 // Pass 3D simple parameters to NETGEN
2805 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
2806 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
2808 if ( double vol = simple3d->GetMaxElementVolume() ) {
2810 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
2811 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2814 // length from faces
2815 mparams.maxh = _ngMesh->AverageH();
2817 _ngMesh->SetGlobalH (mparams.maxh);
2818 mparams.grading = 0.4;
2820 _ngMesh->CalcLocalH(mparams.grading);
2822 _ngMesh->CalcLocalH();
2825 // Care of vertices internal in solids and internal faces (issue 0020676)
2826 if ( internals.hasInternalVertexInSolid() || internals.hasInternalFaces() )
2828 // store computed faces in SMESH in order not to create SMESH
2829 // faces for ng faces added here
2830 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2831 // add ng faces to solids with internal vertices
2832 AddIntVerticesInSolids( occgeo, *_ngMesh, nodeVec, internals );
2833 // duplicate mesh faces on internal faces
2834 FixIntFaces( occgeo, *_ngMesh, internals );
2835 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2837 // Let netgen compute 3D mesh
2838 startWith = endWith = netgen::MESHCONST_MESHVOLUME;
2843 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2845 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2847 if(netgen::multithread.terminate)
2850 if ( comment.empty() ) // do not overwrite a previos error
2851 comment << text(err);
2853 catch (Standard_Failure& ex)
2855 if ( comment.empty() ) // do not overwrite a previos error
2856 comment << text(ex);
2859 catch (netgen::NgException exc)
2861 if ( comment.empty() ) // do not overwrite a previos error
2862 comment << text(exc);
2865 _ticTime = ( doneTime += voluMeshingTime ) / _totalTime / _progressTic;
2867 // Let netgen optimize 3D mesh
2868 if ( !err && _optimize )
2870 startWith = endWith = netgen::MESHCONST_OPTVOLUME;
2875 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2877 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2879 if(netgen::multithread.terminate)
2882 if ( comment.empty() ) // do not overwrite a previos error
2883 comment << text(err);
2885 catch (Standard_Failure& ex)
2887 if ( comment.empty() ) // do not overwrite a previos error
2888 comment << text(ex);
2890 catch (netgen::NgException exc)
2892 if ( comment.empty() ) // do not overwrite a previos error
2893 comment << text(exc);
2897 if (!err && mparams.secondorder > 0)
2902 if ( !meshedSM[ MeshDim_1D ].empty() )
2904 // remove segments not attached to geometry (IPAL0052479)
2905 for (int i = 1; i <= _ngMesh->GetNSeg(); ++i)
2907 const netgen::Segment & seg = _ngMesh->LineSegment (i);
2908 if ( seg.epgeominfo[ 0 ].edgenr == 0 )
2909 _ngMesh->DeleteSegment( i );
2911 _ngMesh->Compress();
2913 // convert to quadratic
2914 netgen::OCCRefinementSurfaces ref (occgeo);
2915 ref.MakeSecondOrder (*_ngMesh);
2917 // care of elements already loaded to SMESH
2918 // if ( initState._nbSegments > 0 )
2919 // makeQuadratic( occgeo.emap, _mesh );
2920 // if ( initState._nbFaces > 0 )
2921 // makeQuadratic( occgeo.fmap, _mesh );
2923 catch (Standard_Failure& ex)
2925 if ( comment.empty() ) // do not overwrite a previos error
2926 comment << "Exception in netgen at passing to 2nd order ";
2928 catch (netgen::NgException exc)
2930 if ( comment.empty() ) // do not overwrite a previos error
2931 comment << exc.What();
2936 _ticTime = 0.98 / _progressTic;
2938 int nbNod = _ngMesh->GetNP();
2939 int nbSeg = _ngMesh->GetNSeg();
2940 int nbFac = _ngMesh->GetNSE();
2941 int nbVol = _ngMesh->GetNE();
2942 bool isOK = ( !err && (_isVolume ? (nbVol > 0) : (nbFac > 0)) );
2944 MESSAGE((err ? "Mesh Generation failure" : "End of Mesh Generation") <<
2945 ", nb nodes: " << nbNod <<
2946 ", nb segments: " << nbSeg <<
2947 ", nb faces: " << nbFac <<
2948 ", nb volumes: " << nbVol);
2950 // Feed back the SMESHDS with the generated Nodes and Elements
2951 if ( true /*isOK*/ ) // get whatever built
2953 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2955 if ( quadHelper.GetIsQuadratic() ) // remove free nodes
2956 for ( size_t i = 0; i < nodeVec.size(); ++i )
2957 if ( nodeVec[i] && nodeVec[i]->NbInverseElements() == 0 )
2958 _mesh->GetMeshDS()->RemoveFreeNode( nodeVec[i], 0, /*fromGroups=*/false );
2960 SMESH_ComputeErrorPtr readErr = ReadErrors(nodeVec);
2961 if ( readErr && !readErr->myBadElements.empty() )
2964 if ( !comment.empty() && !readErr->myComment.empty() ) comment += "\n";
2965 comment += readErr->myComment;
2967 if ( error->IsOK() && ( !isOK || comment.size() > 0 ))
2968 error->myName = COMPERR_ALGO_FAILED;
2969 if ( !comment.empty() )
2970 error->myComment = comment;
2972 // SetIsAlwaysComputed( true ) to empty sub-meshes, which
2973 // appear if the geometry contains coincident sub-shape due
2974 // to bool merge_solids = 1; in netgen/libsrc/occ/occgenmesh.cpp
2975 const int nbMaps = 2;
2976 const TopTools_IndexedMapOfShape* geoMaps[nbMaps] =
2977 { & occgeo.vmap, & occgeo.emap/*, & occgeo.fmap*/ };
2978 for ( int iMap = 0; iMap < nbMaps; ++iMap )
2979 for (int i = 1; i <= geoMaps[iMap]->Extent(); i++)
2980 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( geoMaps[iMap]->FindKey(i)))
2981 if ( !sm->IsMeshComputed() )
2982 sm->SetIsAlwaysComputed( true );
2984 // set bad compute error to subshapes of all failed sub-shapes
2985 if ( !error->IsOK() )
2987 bool pb2D = false, pb3D = false;
2988 for (int i = 1; i <= occgeo.fmap.Extent(); i++) {
2989 int status = occgeo.facemeshstatus[i-1];
2990 if (status == 1 ) continue;
2991 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.fmap( i ))) {
2992 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
2993 if ( !smError || smError->IsOK() ) {
2995 smError.reset( new SMESH_ComputeError( *error ));
2997 smError.reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED, "Ignored" ));
2998 if ( SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
2999 smError->myName = COMPERR_WARNING;
3001 pb2D = pb2D || smError->IsKO();
3004 if ( !pb2D ) // all faces are OK
3005 for (int i = 1; i <= occgeo.somap.Extent(); i++)
3006 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.somap( i )))
3008 bool smComputed = nbVol && !sm->IsEmpty();
3009 if ( smComputed && internals.hasInternalVertexInSolid( sm->GetId() ))
3011 int nbIntV = internals.getSolidsWithVertices().find( sm->GetId() )->second.size();
3012 SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
3013 smComputed = ( smDS->NbElements() > 0 || smDS->NbNodes() > nbIntV );
3015 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
3016 if ( !smComputed && ( !smError || smError->IsOK() ))
3018 smError.reset( new SMESH_ComputeError( *error ));
3019 if ( nbVol && SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
3021 smError->myName = COMPERR_WARNING;
3023 else if ( !smError->myBadElements.empty() ) // bad surface mesh
3025 if ( !hasBadElemOnSolid( smError->myBadElements, sm ))
3029 pb3D = pb3D || ( smError && smError->IsKO() );
3031 if ( !pb2D && !pb3D )
3032 err = 0; // no fatal errors, only warnings
3035 ngLib._isComputeOk = !err;
3040 //=============================================================================
3044 //=============================================================================
3045 bool NETGENPlugin_Mesher::Evaluate(MapShapeNbElems& aResMap)
3047 netgen::MeshingParameters& mparams = netgen::mparam;
3050 // -------------------------
3051 // Prepare OCC geometry
3052 // -------------------------
3053 netgen::OCCGeometry occgeo;
3054 list< SMESH_subMesh* > meshedSM[4]; // for 0-3 dimensions
3055 NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
3056 PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
3058 bool tooManyElems = false;
3059 const int hugeNb = std::numeric_limits<int>::max() / 100;
3064 // pass 1D simple parameters to NETGEN
3067 // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
3068 mparams.uselocalh = false;
3069 mparams.grading = 0.8; // not limitited size growth
3071 if ( _simpleHyp->GetNumberOfSegments() )
3073 mparams.maxh = occgeo.boundingbox.Diam();
3076 mparams.maxh = _simpleHyp->GetLocalLength();
3079 if ( mparams.maxh == 0.0 )
3080 mparams.maxh = occgeo.boundingbox.Diam();
3081 if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
3082 mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
3084 // let netgen create _ngMesh and calculate element size on not meshed shapes
3085 NETGENPlugin_NetgenLibWrapper ngLib;
3086 netgen::Mesh *ngMesh = NULL;
3090 int startWith = netgen::MESHCONST_ANALYSE;
3091 int endWith = netgen::MESHCONST_MESHEDGES;
3093 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, mparams, startWith, endWith);
3095 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, startWith, endWith, optstr);
3098 if(netgen::multithread.terminate)
3101 ngLib.setMesh(( Ng_Mesh*) ngMesh );
3103 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( _shape ))
3104 sm->GetComputeError().reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED ));
3109 // Pass 1D simple parameters to NETGEN
3110 // --------------------------------
3111 int nbSeg = _simpleHyp->GetNumberOfSegments();
3112 double segSize = _simpleHyp->GetLocalLength();
3113 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
3115 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
3117 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
3118 setLocalSize( e, segSize, *ngMesh );
3121 else // if ( ! _simpleHyp )
3123 // Local size on vertices and edges
3124 // --------------------------------
3125 for(std::map<int,double>::const_iterator it=EdgeId2LocalSize.begin(); it!=EdgeId2LocalSize.end(); it++)
3127 int key = (*it).first;
3128 double hi = (*it).second;
3129 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3130 const TopoDS_Edge& e = TopoDS::Edge(shape);
3131 setLocalSize( e, hi, *ngMesh );
3133 for(std::map<int,double>::const_iterator it=VertexId2LocalSize.begin(); it!=VertexId2LocalSize.end(); it++)
3135 int key = (*it).first;
3136 double hi = (*it).second;
3137 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3138 const TopoDS_Vertex& v = TopoDS::Vertex(shape);
3139 gp_Pnt p = BRep_Tool::Pnt(v);
3140 NETGENPlugin_Mesher::RestrictLocalSize( *ngMesh, p.XYZ(), hi );
3142 for(map<int,double>::const_iterator it=FaceId2LocalSize.begin();
3143 it!=FaceId2LocalSize.end(); it++)
3145 int key = (*it).first;
3146 double val = (*it).second;
3147 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3148 int faceNgID = occgeo.fmap.FindIndex(shape);
3149 occgeo.SetFaceMaxH(faceNgID, val);
3150 for ( TopExp_Explorer edgeExp( shape, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
3151 setLocalSize( TopoDS::Edge( edgeExp.Current() ), val, *ngMesh );
3154 // calculate total nb of segments and length of edges
3155 double fullLen = 0.0;
3157 int entity = mparams.secondorder > 0 ? SMDSEntity_Quad_Edge : SMDSEntity_Edge;
3158 TopTools_DataMapOfShapeInteger Edge2NbSeg;
3159 for (TopExp_Explorer exp(_shape, TopAbs_EDGE); exp.More(); exp.Next())
3161 TopoDS_Edge E = TopoDS::Edge( exp.Current() );
3162 if( !Edge2NbSeg.Bind(E,0) )
3165 double aLen = SMESH_Algo::EdgeLength(E);
3168 vector<int>& aVec = aResMap[_mesh->GetSubMesh(E)];
3170 aVec.resize( SMDSEntity_Last, 0);
3172 fullNbSeg += aVec[ entity ];
3175 // store nb of segments computed by Netgen
3176 NCollection_Map<Link> linkMap;
3177 for (int i = 1; i <= ngMesh->GetNSeg(); ++i )
3179 const netgen::Segment& seg = ngMesh->LineSegment(i);
3180 Link link(seg[0], seg[1]);
3181 if ( !linkMap.Add( link )) continue;
3182 int aGeomEdgeInd = seg.epgeominfo[0].edgenr;
3183 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
3185 vector<int>& aVec = aResMap[_mesh->GetSubMesh(occgeo.emap(aGeomEdgeInd))];
3189 // store nb of nodes on edges computed by Netgen
3190 TopTools_DataMapIteratorOfDataMapOfShapeInteger Edge2NbSegIt(Edge2NbSeg);
3191 for (; Edge2NbSegIt.More(); Edge2NbSegIt.Next())
3193 vector<int>& aVec = aResMap[_mesh->GetSubMesh(Edge2NbSegIt.Key())];
3194 if ( aVec[ entity ] > 1 && aVec[ SMDSEntity_Node ] == 0 )
3195 aVec[SMDSEntity_Node] = mparams.secondorder > 0 ? 2*aVec[ entity ]-1 : aVec[ entity ]-1;
3197 fullNbSeg += aVec[ entity ];
3198 Edge2NbSeg( Edge2NbSegIt.Key() ) = aVec[ entity ];
3200 if ( fullNbSeg == 0 )
3207 if ( double area = _simpleHyp->GetMaxElementArea() ) {
3209 mparams.maxh = sqrt(2. * area/sqrt(3.0));
3210 mparams.grading = 0.4; // moderate size growth
3213 // length from edges
3214 mparams.maxh = fullLen/fullNbSeg;
3215 mparams.grading = 0.2; // slow size growth
3218 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3219 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3221 for (TopExp_Explorer exp(_shape, TopAbs_FACE); exp.More(); exp.Next())
3223 TopoDS_Face F = TopoDS::Face( exp.Current() );
3224 SMESH_subMesh *sm = _mesh->GetSubMesh(F);
3226 BRepGProp::SurfaceProperties(F,G);
3227 double anArea = G.Mass();
3228 tooManyElems = tooManyElems || ( anArea/hugeNb > mparams.maxh*mparams.maxh );
3230 if ( !tooManyElems )
3232 TopTools_MapOfShape egdes;
3233 for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next())
3234 if ( egdes.Add( exp1.Current() ))
3235 nb1d += Edge2NbSeg.Find(exp1.Current());
3237 int nbFaces = tooManyElems ? hugeNb : int( 4*anArea / (mparams.maxh*mparams.maxh*sqrt(3.)));
3238 int nbNodes = tooManyElems ? hugeNb : (( nbFaces*3 - (nb1d-1)*2 ) / 6 + 1 );
3240 vector<int> aVec(SMDSEntity_Last, 0);
3241 if( mparams.secondorder > 0 ) {
3242 int nb1d_in = (nbFaces*3 - nb1d) / 2;
3243 aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
3244 aVec[SMDSEntity_Quad_Triangle] = nbFaces;
3247 aVec[SMDSEntity_Node] = Max ( nbNodes, 0 );
3248 aVec[SMDSEntity_Triangle] = nbFaces;
3250 aResMap[sm].swap(aVec);
3257 // pass 3D simple parameters to NETGEN
3258 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
3259 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
3261 if ( double vol = simple3d->GetMaxElementVolume() ) {
3263 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
3264 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3267 // using previous length from faces
3269 mparams.grading = 0.4;
3270 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3273 BRepGProp::VolumeProperties(_shape,G);
3274 double aVolume = G.Mass();
3275 double tetrVol = 0.1179*mparams.maxh*mparams.maxh*mparams.maxh;
3276 tooManyElems = tooManyElems || ( aVolume/hugeNb > tetrVol );
3277 int nbVols = tooManyElems ? hugeNb : int(aVolume/tetrVol);
3278 int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
3279 vector<int> aVec(SMDSEntity_Last, 0 );
3280 if ( tooManyElems ) // avoid FPE
3282 aVec[SMDSEntity_Node] = hugeNb;
3283 aVec[ mparams.secondorder > 0 ? SMDSEntity_Quad_Tetra : SMDSEntity_Tetra] = hugeNb;
3287 if( mparams.secondorder > 0 ) {
3288 aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
3289 aVec[SMDSEntity_Quad_Tetra] = nbVols;
3292 aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
3293 aVec[SMDSEntity_Tetra] = nbVols;
3296 SMESH_subMesh *sm = _mesh->GetSubMesh(_shape);
3297 aResMap[sm].swap(aVec);
3303 double NETGENPlugin_Mesher::GetProgress(const SMESH_Algo* holder,
3304 const int * algoProgressTic,
3305 const double * algoProgress) const
3307 ((int&) _progressTic ) = *algoProgressTic + 1;
3309 if ( !_occgeom ) return 0;
3311 double progress = -1;
3314 if ( _ticTime < 0 && netgen::multithread.task[0] == 'O'/*Optimizing surface*/ )
3316 ((double&) _ticTime ) = edgeFaceMeshingTime / _totalTime / _progressTic;
3318 else if ( !_optimize /*&& _occgeom->fmap.Extent() > 1*/ )
3320 int doneShapeIndex = -1;
3321 while ( doneShapeIndex+1 < _occgeom->facemeshstatus.Size() &&
3322 _occgeom->facemeshstatus[ doneShapeIndex+1 ])
3324 if ( doneShapeIndex+1 != _curShapeIndex )
3326 ((int&) _curShapeIndex) = doneShapeIndex+1;
3327 double doneShapeRate = _curShapeIndex / double( _occgeom->fmap.Extent() );
3328 double doneTime = edgeMeshingTime + doneShapeRate * faceMeshingTime;
3329 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3330 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3331 // << " " << doneTime / _totalTime / _progressTic << endl;
3335 else if ( !_optimize && _occgeom->somap.Extent() > 1 )
3337 int curShapeIndex = _curShapeIndex;
3338 if ( _ngMesh->GetNE() > 0 )
3340 netgen::Element el = (*_ngMesh)[netgen::ElementIndex( _ngMesh->GetNE()-1 )];
3341 curShapeIndex = el.GetIndex();
3343 if ( curShapeIndex != _curShapeIndex )
3345 ((int&) _curShapeIndex) = curShapeIndex;
3346 double doneShapeRate = _curShapeIndex / double( _occgeom->somap.Extent() );
3347 double doneTime = edgeFaceMeshingTime + doneShapeRate * voluMeshingTime;
3348 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3349 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3350 // << " " << doneTime / _totalTime / _progressTic << endl;
3354 progress = Max( *algoProgressTic * _ticTime, *algoProgress );
3357 ((int&) *algoProgressTic )++;
3358 ((double&) *algoProgress) = progress;
3360 //cout << progress << " " << *algoProgressTic << " " << netgen::multithread.task << " "<< _ticTime << endl;
3362 return Min( progress, 0.99 );
3365 //================================================================================
3367 * \brief Remove "test.out" and "problemfaces" files in current directory
3369 //================================================================================
3371 void NETGENPlugin_Mesher::RemoveTmpFiles()
3373 bool rm = SMESH_File("test.out").remove() ;
3375 if (rm && netgen::testout)
3377 delete netgen::testout;
3378 netgen::testout = 0;
3381 SMESH_File("problemfaces").remove();
3382 SMESH_File("occmesh.rep").remove();
3385 //================================================================================
3387 * \brief Read mesh entities preventing successful computation from "test.out" file
3389 //================================================================================
3391 SMESH_ComputeErrorPtr
3392 NETGENPlugin_Mesher::ReadErrors(const vector<const SMDS_MeshNode* >& nodeVec)
3394 SMESH_ComputeErrorPtr err = SMESH_ComputeError::New
3395 (COMPERR_BAD_INPUT_MESH, "Some edges multiple times in surface mesh");
3396 SMESH_File file("test.out");
3398 vector<int> three1(3), three2(3);
3399 const char* badEdgeStr = " multiple times in surface mesh";
3400 const int badEdgeStrLen = strlen( badEdgeStr );
3401 const int nbNodes = nodeVec.size();
3403 while( !file.eof() )
3405 if ( strncmp( file, "Edge ", 5 ) == 0 &&
3406 file.getInts( two ) &&
3407 strncmp( file, badEdgeStr, badEdgeStrLen ) == 0 &&
3408 two[0] < nbNodes && two[1] < nbNodes )
3410 err->myBadElements.push_back( new SMDS_LinearEdge( nodeVec[ two[0]], nodeVec[ two[1]] ));
3411 file += badEdgeStrLen;
3413 else if ( strncmp( file, "Intersecting: ", 14 ) == 0 )
3416 // openelement 18 with open element 126
3420 const char* pos = file;
3421 bool ok = ( strncmp( file, "openelement ", 12 ) == 0 );
3422 ok = ok && file.getInts( two );
3423 ok = ok && file.getInts( three1 );
3424 ok = ok && file.getInts( three2 );
3425 for ( int i = 0; ok && i < 3; ++i )
3426 ok = ( three1[i] < nbNodes && nodeVec[ three1[i]]);
3427 for ( int i = 0; ok && i < 3; ++i )
3428 ok = ( three2[i] < nbNodes && nodeVec[ three2[i]]);
3431 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three1[0]],
3432 nodeVec[ three1[1]],
3433 nodeVec[ three1[2]]));
3434 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three2[0]],
3435 nodeVec[ three2[1]],
3436 nodeVec[ three2[2]]));
3437 err->myComment = "Intersecting triangles";
3451 size_t nbBadElems = err->myBadElements.size();
3458 //================================================================================
3460 * \brief Write a python script creating an equivalent SALOME mesh.
3461 * This is useful to see what mesh is passed as input for the next step of mesh
3462 * generation (of mesh of higher dimension)
3464 //================================================================================
3466 void NETGENPlugin_Mesher::toPython( const netgen::Mesh* ngMesh,
3467 const std::string& pyFile)
3469 ofstream outfile(pyFile.c_str(), ios::out);
3470 if ( !outfile ) return;
3472 outfile << "import SMESH" << endl
3473 << "from salome.smesh import smeshBuilder" << endl
3474 << "smesh = smeshBuilder.New(salome.myStudy)" << endl
3475 << "mesh = smesh.Mesh()" << endl << endl;
3477 using namespace netgen;
3479 for (pi = PointIndex::BASE;
3480 pi < ngMesh->GetNP()+PointIndex::BASE; pi++)
3482 outfile << "mesh.AddNode( ";
3483 outfile << (*ngMesh)[pi](0) << ", ";
3484 outfile << (*ngMesh)[pi](1) << ", ";
3485 outfile << (*ngMesh)[pi](2) << ") ## "<< pi << endl;
3488 int nbDom = ngMesh->GetNDomains();
3489 for ( int i = 0; i < nbDom; ++i )
3490 outfile<< "grp" << i+1 << " = mesh.CreateEmptyGroup( SMESH.FACE, 'domain"<< i+1 << "')"<< endl;
3492 SurfaceElementIndex sei;
3493 for (sei = 0; sei < ngMesh->GetNSE(); sei++)
3495 outfile << "mesh.AddFace([ ";
3496 Element2d sel = (*ngMesh)[sei];
3497 for (int j = 0; j < sel.GetNP(); j++)
3498 outfile << sel[j] << ( j+1 < sel.GetNP() ? ", " : " ])");
3499 if ( sel.IsDeleted() ) outfile << " ## IsDeleted ";
3502 if ((*ngMesh)[sei].GetIndex())
3504 if ( int dom1 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainIn())
3505 outfile << "grp"<< dom1 <<".Add([ " << (int)sei+1 << " ])" << endl;
3506 if ( int dom2 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainOut())
3507 outfile << "grp"<< dom2 <<".Add([ " << (int)sei+1 << " ])" << endl;
3511 for (ElementIndex ei = 0; ei < ngMesh->GetNE(); ei++)
3513 Element el = (*ngMesh)[ei];
3514 outfile << "mesh.AddVolume([ ";
3515 for (int j = 0; j < el.GetNP(); j++)
3516 outfile << el[j] << ( j+1 < el.GetNP() ? ", " : " ])");
3520 for (int i = 1; i <= ngMesh->GetNSeg(); i++)
3522 const Segment & seg = ngMesh->LineSegment (i);
3523 outfile << "mesh.AddEdge([ "
3525 << seg[1] << " ])" << endl;
3527 cout << "Write " << pyFile << endl;
3530 //================================================================================
3532 * \brief Constructor of NETGENPlugin_ngMeshInfo
3534 //================================================================================
3536 NETGENPlugin_ngMeshInfo::NETGENPlugin_ngMeshInfo( netgen::Mesh* ngMesh):
3541 _nbNodes = ngMesh->GetNP();
3542 _nbSegments = ngMesh->GetNSeg();
3543 _nbFaces = ngMesh->GetNSE();
3544 _nbVolumes = ngMesh->GetNE();
3548 _nbNodes = _nbSegments = _nbFaces = _nbVolumes = 0;
3552 //================================================================================
3554 * \brief Copy LocalH member from one netgen mesh to another
3556 //================================================================================
3558 void NETGENPlugin_ngMeshInfo::transferLocalH( netgen::Mesh* fromMesh,
3559 netgen::Mesh* toMesh )
3561 if ( !fromMesh->LocalHFunctionGenerated() ) return;
3562 if ( !toMesh->LocalHFunctionGenerated() )
3564 toMesh->CalcLocalH(netgen::mparam.grading);
3566 toMesh->CalcLocalH();
3569 const size_t size = sizeof( netgen::LocalH );
3570 _copyOfLocalH = new char[ size ];
3571 memcpy( (void*)_copyOfLocalH, (void*)&toMesh->LocalHFunction(), size );
3572 memcpy( (void*)&toMesh->LocalHFunction(), (void*)&fromMesh->LocalHFunction(), size );
3575 //================================================================================
3577 * \brief Restore LocalH member of a netgen mesh
3579 //================================================================================
3581 void NETGENPlugin_ngMeshInfo::restoreLocalH( netgen::Mesh* toMesh )
3583 if ( _copyOfLocalH )
3585 const size_t size = sizeof( netgen::LocalH );
3586 memcpy( (void*)&toMesh->LocalHFunction(), (void*)_copyOfLocalH, size );
3587 delete [] _copyOfLocalH;
3592 //================================================================================
3594 * \brief Find "internal" sub-shapes
3596 //================================================================================
3598 NETGENPlugin_Internals::NETGENPlugin_Internals( SMESH_Mesh& mesh,
3599 const TopoDS_Shape& shape,
3601 : _mesh( mesh ), _is3D( is3D )
3603 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3605 TopExp_Explorer f,e;
3606 for ( f.Init( shape, TopAbs_FACE ); f.More(); f.Next() )
3608 int faceID = meshDS->ShapeToIndex( f.Current() );
3610 // find not computed internal edges
3612 for ( e.Init( f.Current().Oriented(TopAbs_FORWARD), TopAbs_EDGE ); e.More(); e.Next() )
3613 if ( e.Current().Orientation() == TopAbs_INTERNAL )
3615 SMESH_subMesh* eSM = mesh.GetSubMesh( e.Current() );
3616 if ( eSM->IsEmpty() )
3618 _e2face.insert( make_pair( eSM->GetId(), faceID ));
3619 for ( TopoDS_Iterator v(e.Current()); v.More(); v.Next() )
3620 _e2face.insert( make_pair( meshDS->ShapeToIndex( v.Value() ), faceID ));
3624 // find internal vertices in a face
3625 set<int> intVV; // issue 0020850 where same vertex is twice in a face
3626 for ( TopoDS_Iterator fSub( f.Current() ); fSub.More(); fSub.Next())
3627 if ( fSub.Value().ShapeType() == TopAbs_VERTEX )
3629 int vID = meshDS->ShapeToIndex( fSub.Value() );
3630 if ( intVV.insert( vID ).second )
3631 _f2v[ faceID ].push_back( vID );
3636 // find internal faces and their subshapes where nodes are to be doubled
3637 // to make a crack with non-sewed borders
3639 if ( f.Current().Orientation() == TopAbs_INTERNAL )
3641 _intShapes.insert( meshDS->ShapeToIndex( f.Current() ));
3644 list< TopoDS_Shape > edges;
3645 for ( e.Init( f.Current(), TopAbs_EDGE ); e.More(); e.Next())
3646 if ( SMESH_MesherHelper::NbAncestors( e.Current(), mesh, TopAbs_FACE ) > 1 )
3648 _intShapes.insert( meshDS->ShapeToIndex( e.Current() ));
3649 edges.push_back( e.Current() );
3650 // find border faces
3651 PShapeIteratorPtr fIt =
3652 SMESH_MesherHelper::GetAncestors( edges.back(),mesh,TopAbs_FACE );
3653 while ( const TopoDS_Shape* pFace = fIt->next() )
3654 if ( !pFace->IsSame( f.Current() ))
3655 _borderFaces.insert( meshDS->ShapeToIndex( *pFace ));
3658 // we consider vertex internal if it is shared by more than one internal edge
3659 list< TopoDS_Shape >::iterator edge = edges.begin();
3660 for ( ; edge != edges.end(); ++edge )
3661 for ( TopoDS_Iterator v( *edge ); v.More(); v.Next() )
3663 set<int> internalEdges;
3664 PShapeIteratorPtr eIt =
3665 SMESH_MesherHelper::GetAncestors( v.Value(),mesh,TopAbs_EDGE );
3666 while ( const TopoDS_Shape* pEdge = eIt->next() )
3668 int edgeID = meshDS->ShapeToIndex( *pEdge );
3669 if ( isInternalShape( edgeID ))
3670 internalEdges.insert( edgeID );
3672 if ( internalEdges.size() > 1 )
3673 _intShapes.insert( meshDS->ShapeToIndex( v.Value() ));
3677 } // loop on geom faces
3679 // find vertices internal in solids
3682 for ( TopExp_Explorer so(shape, TopAbs_SOLID); so.More(); so.Next())
3684 int soID = meshDS->ShapeToIndex( so.Current() );
3685 for ( TopoDS_Iterator soSub( so.Current() ); soSub.More(); soSub.Next())
3686 if ( soSub.Value().ShapeType() == TopAbs_VERTEX )
3687 _s2v[ soID ].push_back( meshDS->ShapeToIndex( soSub.Value() ));
3692 //================================================================================
3694 * \brief Find mesh faces on non-internal geom faces sharing internal edge
3695 * some nodes of which are to be doubled to make the second border of the "crack"
3697 //================================================================================
3699 void NETGENPlugin_Internals::findBorderElements( TIDSortedElemSet & borderElems )
3701 if ( _intShapes.empty() ) return;
3703 SMESH_Mesh& mesh = const_cast<SMESH_Mesh&>(_mesh);
3704 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3706 // loop on internal geom edges
3707 set<int>::const_iterator intShapeId = _intShapes.begin();
3708 for ( ; intShapeId != _intShapes.end(); ++intShapeId )
3710 const TopoDS_Shape& s = meshDS->IndexToShape( *intShapeId );
3711 if ( s.ShapeType() != TopAbs_EDGE ) continue;
3713 // get internal and non-internal geom faces sharing the internal edge <s>
3715 set<int>::iterator bordFace = _borderFaces.end();
3716 PShapeIteratorPtr faces = SMESH_MesherHelper::GetAncestors( s, _mesh, TopAbs_FACE );
3717 while ( const TopoDS_Shape* pFace = faces->next() )
3719 int faceID = meshDS->ShapeToIndex( *pFace );
3720 if ( isInternalShape( faceID ))
3723 bordFace = _borderFaces.insert( faceID ).first;
3725 if ( bordFace == _borderFaces.end() || !intFace ) continue;
3727 // get all links of mesh faces on internal geom face sharing nodes on edge <s>
3728 set< SMESH_OrientedLink > links; //!< links of faces on internal geom face
3729 list<const SMDS_MeshElement*> suspectFaces[2]; //!< mesh faces on border geom faces
3730 int nbSuspectFaces = 0;
3731 SMESHDS_SubMesh* intFaceSM = meshDS->MeshElements( intFace );
3732 if ( !intFaceSM || intFaceSM->NbElements() == 0 ) continue;
3733 SMESH_subMeshIteratorPtr smIt = mesh.GetSubMesh( s )->getDependsOnIterator(true,true);
3734 while ( smIt->more() )
3736 SMESHDS_SubMesh* sm = smIt->next()->GetSubMeshDS();
3737 if ( !sm ) continue;
3738 SMDS_NodeIteratorPtr nIt = sm->GetNodes();
3739 while ( nIt->more() )
3741 const SMDS_MeshNode* nOnEdge = nIt->next();
3742 SMDS_ElemIteratorPtr fIt = nOnEdge->GetInverseElementIterator(SMDSAbs_Face);
3743 while ( fIt->more() )
3745 const SMDS_MeshElement* f = fIt->next();
3746 const int nbNodes = f->NbCornerNodes();
3747 if ( intFaceSM->Contains( f ))
3749 for ( int i = 0; i < nbNodes; ++i )
3750 links.insert( SMESH_OrientedLink( f->GetNode(i), f->GetNode((i+1)%nbNodes)));
3755 for ( int i = 0; i < nbNodes; ++i )
3756 nbDblNodes += isInternalShape( f->GetNode(i)->getshapeId() );
3758 suspectFaces[ nbDblNodes < 2 ].push_back( f );
3764 // suspectFaces[0] having link with same orientation as mesh faces on
3765 // the internal geom face are <borderElems>. suspectFaces[1] have
3766 // only one node on edge <s>, we decide on them later (at the 2nd loop)
3767 // by links of <borderElems> found at the 1st and 2nd loops
3768 set< SMESH_OrientedLink > borderLinks;
3769 for ( int isPostponed = 0; isPostponed < 2; ++isPostponed )
3771 list<const SMDS_MeshElement*>::iterator fIt = suspectFaces[isPostponed].begin();
3772 for ( int nbF = 0; fIt != suspectFaces[isPostponed].end(); ++fIt, ++nbF )
3774 const SMDS_MeshElement* f = *fIt;
3775 bool isBorder = false, linkFound = false, borderLinkFound = false;
3776 list< SMESH_OrientedLink > faceLinks;
3777 int nbNodes = f->NbCornerNodes();
3778 for ( int i = 0; i < nbNodes; ++i )
3780 SMESH_OrientedLink link( f->GetNode(i), f->GetNode((i+1)%nbNodes));
3781 faceLinks.push_back( link );
3784 set< SMESH_OrientedLink >::iterator foundLink = links.find( link );
3785 if ( foundLink != links.end() )
3788 isBorder = ( foundLink->_reversed == link._reversed );
3789 if ( !isBorder && !isPostponed ) break;
3790 faceLinks.pop_back();
3792 else if ( isPostponed && !borderLinkFound )
3794 foundLink = borderLinks.find( link );
3795 if ( foundLink != borderLinks.end() )
3797 borderLinkFound = true;
3798 isBorder = ( foundLink->_reversed != link._reversed );
3805 borderElems.insert( f );
3806 borderLinks.insert( faceLinks.begin(), faceLinks.end() );
3808 else if ( !linkFound && !borderLinkFound )
3810 suspectFaces[1].push_back( f );
3811 if ( nbF > 2 * nbSuspectFaces )
3812 break; // dead loop protection
3819 //================================================================================
3821 * \brief put internal shapes in maps and fill in submeshes to precompute
3823 //================================================================================
3825 void NETGENPlugin_Internals::getInternalEdges( TopTools_IndexedMapOfShape& fmap,
3826 TopTools_IndexedMapOfShape& emap,
3827 TopTools_IndexedMapOfShape& vmap,
3828 list< SMESH_subMesh* > smToPrecompute[])
3830 if ( !hasInternalEdges() ) return;
3831 map<int,int>::const_iterator ev_face = _e2face.begin();
3832 for ( ; ev_face != _e2face.end(); ++ev_face )
3834 const TopoDS_Shape& ev = _mesh.GetMeshDS()->IndexToShape( ev_face->first );
3835 const TopoDS_Shape& face = _mesh.GetMeshDS()->IndexToShape( ev_face->second );
3837 ( ev.ShapeType() == TopAbs_EDGE ? emap : vmap ).Add( ev );
3839 //cout<<"INTERNAL EDGE or VERTEX "<<ev_face->first<<" on face "<<ev_face->second<<endl;
3841 smToPrecompute[ MeshDim_1D ].push_back( _mesh.GetSubMeshContaining( ev_face->first ));
3845 //================================================================================
3847 * \brief return shapes and submeshes to be meshed and already meshed boundary submeshes
3849 //================================================================================
3851 void NETGENPlugin_Internals::getInternalFaces( TopTools_IndexedMapOfShape& fmap,
3852 TopTools_IndexedMapOfShape& emap,
3853 list< SMESH_subMesh* >& intFaceSM,
3854 list< SMESH_subMesh* >& boundarySM)
3856 if ( !hasInternalFaces() ) return;
3858 // <fmap> and <emap> are for not yet meshed shapes
3859 // <intFaceSM> is for submeshes of faces
3860 // <boundarySM> is for meshed edges and vertices
3865 set<int> shapeIDs ( _intShapes );
3866 if ( !_borderFaces.empty() )
3867 shapeIDs.insert( _borderFaces.begin(), _borderFaces.end() );
3869 set<int>::const_iterator intS = shapeIDs.begin();
3870 for ( ; intS != shapeIDs.end(); ++intS )
3872 SMESH_subMesh* sm = _mesh.GetSubMeshContaining( *intS );
3874 if ( sm->GetSubShape().ShapeType() != TopAbs_FACE ) continue;
3876 intFaceSM.push_back( sm );
3878 // add submeshes of not computed internal faces
3879 if ( !sm->IsEmpty() ) continue;
3881 SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(true,true);
3882 while ( smIt->more() )
3885 const TopoDS_Shape& s = sm->GetSubShape();
3887 if ( sm->IsEmpty() )
3890 switch ( s.ShapeType() ) {
3891 case TopAbs_FACE: fmap.Add ( s ); break;
3892 case TopAbs_EDGE: emap.Add ( s ); break;
3898 if ( s.ShapeType() != TopAbs_FACE )
3899 boundarySM.push_back( sm );
3905 //================================================================================
3907 * \brief Return true if given shape is to be precomputed in order to be correctly
3908 * added to netgen mesh
3910 //================================================================================
3912 bool NETGENPlugin_Internals::isShapeToPrecompute(const TopoDS_Shape& s)
3914 int shapeID = _mesh.GetMeshDS()->ShapeToIndex( s );
3915 switch ( s.ShapeType() ) {
3916 case TopAbs_FACE : break; //return isInternalShape( shapeID ) || isBorderFace( shapeID );
3917 case TopAbs_EDGE : return isInternalEdge( shapeID );
3918 case TopAbs_VERTEX: break;
3924 //================================================================================
3926 * \brief Return SMESH
3928 //================================================================================
3930 SMESH_Mesh& NETGENPlugin_Internals::getMesh() const
3932 return const_cast<SMESH_Mesh&>( _mesh );
3935 //================================================================================
3937 * \brief Initialize netgen library
3939 //================================================================================
3941 NETGENPlugin_NetgenLibWrapper::NETGENPlugin_NetgenLibWrapper()
3945 _isComputeOk = false;
3947 if ( !getenv( "KEEP_NETGEN_OUTPUT" ))
3949 // redirect all netgen output (mycout,myerr,cout) to _outputFileName
3950 _outputFileName = getOutputFileName();
3951 netgen::mycout = new ofstream ( _outputFileName.c_str() );
3952 netgen::myerr = netgen::mycout;
3953 _coutBuffer = std::cout.rdbuf();
3955 cout << "NOTE: netgen output is redirected to file " << _outputFileName << endl;
3957 std::cout.rdbuf( netgen::mycout->rdbuf() );
3961 _ngMesh = Ng_NewMesh();
3964 //================================================================================
3966 * \brief Finish using netgen library
3968 //================================================================================
3970 NETGENPlugin_NetgenLibWrapper::~NETGENPlugin_NetgenLibWrapper()
3972 Ng_DeleteMesh( _ngMesh );
3974 NETGENPlugin_Mesher::RemoveTmpFiles();
3976 std::cout.rdbuf( _coutBuffer );
3983 //================================================================================
3985 * \brief Set netgen mesh to delete at destruction
3987 //================================================================================
3989 void NETGENPlugin_NetgenLibWrapper::setMesh( Ng_Mesh* mesh )
3992 Ng_DeleteMesh( _ngMesh );
3996 //================================================================================
3998 * \brief Return a unique file name
4000 //================================================================================
4002 std::string NETGENPlugin_NetgenLibWrapper::getOutputFileName()
4004 std::string aTmpDir = SALOMEDS_Tool::GetTmpDir();
4006 TCollection_AsciiString aGenericName = (char*)aTmpDir.c_str();
4007 aGenericName += "NETGEN_";
4009 aGenericName += getpid();
4011 aGenericName += _getpid();
4013 aGenericName += "_";
4014 aGenericName += Abs((Standard_Integer)(long) aGenericName.ToCString());
4015 aGenericName += ".out";
4017 return aGenericName.ToCString();
4020 //================================================================================
4022 * \brief Remove file with netgen output
4024 //================================================================================
4026 void NETGENPlugin_NetgenLibWrapper::removeOutputFile()
4028 if ( !_outputFileName.empty() )
4030 if ( netgen::mycout )
4032 delete netgen::mycout;
4036 string tmpDir = SALOMEDS_Tool::GetDirFromPath ( _outputFileName );
4037 string aFileName = SALOMEDS_Tool::GetNameFromPath( _outputFileName ) + ".out";
4038 SALOMEDS::ListOfFileNames_var aFiles = new SALOMEDS::ListOfFileNames;
4040 aFiles[0] = aFileName.c_str();
4042 SALOMEDS_Tool::RemoveTemporaryFiles( tmpDir.c_str(), aFiles.in(), true );