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 ; 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 ( overrideMinH && netgen::mparam.minh > size )
710 ngMesh.SetMinimalH( size );
711 netgen::mparam.minh = size;
713 netgen::Point3d pi(p.X(), p.Y(), p.Z());
714 ngMesh.RestrictLocalH( pi, size );
717 //================================================================================
719 * \brief fill ngMesh with nodes and elements of computed submeshes
721 //================================================================================
723 bool NETGENPlugin_Mesher::FillNgMesh(netgen::OCCGeometry& occgeom,
724 netgen::Mesh& ngMesh,
725 vector<const SMDS_MeshNode*>& nodeVec,
726 const list< SMESH_subMesh* > & meshedSM,
727 SMESH_MesherHelper* quadHelper,
728 SMESH_ProxyMesh::Ptr proxyMesh)
730 TNode2IdMap nodeNgIdMap;
731 for ( int i = 1; i < nodeVec.size(); ++i )
732 nodeNgIdMap.insert( make_pair( nodeVec[i], i ));
734 TopTools_MapOfShape visitedShapes;
735 map< SMESH_subMesh*, set< int > > visitedEdgeSM2Faces;
736 set< SMESH_subMesh* > computedSM( meshedSM.begin(), meshedSM.end() );
738 SMESH_MesherHelper helper (*_mesh);
740 int faceNgID = ngMesh.GetNFD();
742 list< SMESH_subMesh* >::const_iterator smIt, smEnd = meshedSM.end();
743 for ( smIt = meshedSM.begin(); smIt != smEnd; ++smIt )
745 SMESH_subMesh* sm = *smIt;
746 if ( !visitedShapes.Add( sm->GetSubShape() ))
749 const SMESHDS_SubMesh * smDS = sm->GetSubMeshDS();
750 if ( !smDS ) continue;
752 switch ( sm->GetSubShape().ShapeType() )
754 case TopAbs_EDGE: { // EDGE
755 // ----------------------
756 TopoDS_Edge geomEdge = TopoDS::Edge( sm->GetSubShape() );
757 if ( geomEdge.Orientation() >= TopAbs_INTERNAL )
758 geomEdge.Orientation( TopAbs_FORWARD ); // issue 0020676
760 // Add ng segments for each not meshed FACE the EDGE bounds
761 PShapeIteratorPtr fIt = helper.GetAncestors( geomEdge, *sm->GetFather(), TopAbs_FACE );
762 while ( const TopoDS_Shape * anc = fIt->next() )
764 faceNgID = occgeom.fmap.FindIndex( *anc );
766 continue; // meshed face
768 int faceSMDSId = helper.GetMeshDS()->ShapeToIndex( *anc );
769 if ( visitedEdgeSM2Faces[ sm ].count( faceSMDSId ))
770 continue; // already treated EDGE
772 TopoDS_Face face = TopoDS::Face( occgeom.fmap( faceNgID ));
773 if ( face.Orientation() >= TopAbs_INTERNAL )
774 face.Orientation( TopAbs_FORWARD ); // issue 0020676
776 // get all meshed EDGEs of the FACE connected to geomEdge (issue 0021140)
777 helper.SetSubShape( face );
778 list< TopoDS_Edge > edges = getConnectedEdges( geomEdge, face, computedSM, helper,
779 visitedEdgeSM2Faces );
781 continue; // wrong ancestor?
783 // find out orientation of <edges> within <face>
784 TopoDS_Edge eNotSeam = edges.front();
785 if ( helper.HasSeam() )
787 list< TopoDS_Edge >::iterator eIt = edges.begin();
788 while ( helper.IsRealSeam( *eIt )) ++eIt;
789 if ( eIt != edges.end() )
792 TopAbs_Orientation fOri = helper.GetSubShapeOri( face, eNotSeam );
793 bool isForwad = ( fOri == eNotSeam.Orientation() || fOri >= TopAbs_INTERNAL );
795 // get all nodes from connected <edges>
796 const bool isQuad = smDS->IsQuadratic();
797 StdMeshers_FaceSide fSide( face, edges, _mesh, isForwad, isQuad );
798 const vector<UVPtStruct>& points = fSide.GetUVPtStruct();
799 if ( points.empty() )
800 return false; // invalid node params?
801 int i, nbSeg = fSide.NbSegments();
803 // remember EDGEs of fSide to treat only once
804 for ( int iE = 0; iE < fSide.NbEdges(); ++iE )
805 visitedEdgeSM2Faces[ helper.GetMesh()->GetSubMesh( fSide.Edge(iE )) ].insert(faceSMDSId);
807 double otherSeamParam = 0;
812 int prevNgId = ngNodeId( points[0].node, ngMesh, nodeNgIdMap );
814 for ( i = 0; i < nbSeg; ++i )
816 const UVPtStruct& p1 = points[ i ];
817 const UVPtStruct& p2 = points[ i+1 ];
819 if ( p1.node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX ) //an EDGE begins
822 if ( helper.IsRealSeam( p1.node->getshapeId() ))
824 TopoDS_Edge e = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
825 isSeam = helper.IsRealSeam( e );
828 otherSeamParam = helper.GetOtherParam( helper.GetPeriodicIndex() & 1 ? p2.u : p2.v );
835 seg[1] = prevNgId = ngNodeId( p2.node, ngMesh, nodeNgIdMap );
836 // node param on curve
837 seg.epgeominfo[ 0 ].dist = p1.param;
838 seg.epgeominfo[ 1 ].dist = p2.param;
840 seg.epgeominfo[ 0 ].u = p1.u;
841 seg.epgeominfo[ 0 ].v = p1.v;
842 seg.epgeominfo[ 1 ].u = p2.u;
843 seg.epgeominfo[ 1 ].v = p2.v;
845 //geomEdge = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
846 //seg.epgeominfo[ 0 ].edgenr = seg.epgeominfo[ 1 ].edgenr = occgeom.emap.FindIndex( geomEdge );
848 //seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
849 seg.si = faceNgID; // = geom.fmap.FindIndex (face);
850 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
851 ngMesh.AddSegment (seg);
853 SMESH_TNodeXYZ np1( p1.node ), np2( p2.node );
854 RestrictLocalSize( ngMesh, 0.5*(np1+np2), (np1-np2).Modulus() );
857 cout << "Segment: " << seg.edgenr << " on SMESH face " << helper.GetMeshDS()->ShapeToIndex( face ) << endl
858 << "\tface index: " << seg.si << endl
859 << "\tp1: " << seg[0] << endl
860 << "\tp2: " << seg[1] << endl
861 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
862 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
863 //<< "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
864 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
865 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl;
866 //<< "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
870 if ( helper.GetPeriodicIndex() && 1 ) {
871 seg.epgeominfo[ 0 ].u = otherSeamParam;
872 seg.epgeominfo[ 1 ].u = otherSeamParam;
873 swap (seg.epgeominfo[0].v, seg.epgeominfo[1].v);
875 seg.epgeominfo[ 0 ].v = otherSeamParam;
876 seg.epgeominfo[ 1 ].v = otherSeamParam;
877 swap (seg.epgeominfo[0].u, seg.epgeominfo[1].u);
879 swap (seg[0], seg[1]);
880 swap (seg.epgeominfo[0].dist, seg.epgeominfo[1].dist);
881 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
882 ngMesh.AddSegment (seg);
884 cout << "Segment: " << seg.edgenr << endl
885 << "\t is SEAM (reverse) of the previous. "
886 << " Other " << (helper.GetPeriodicIndex() && 1 ? "U" : "V")
887 << " = " << otherSeamParam << endl;
890 else if ( fOri == TopAbs_INTERNAL )
892 swap (seg[0], seg[1]);
893 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
894 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
895 ngMesh.AddSegment (seg);
897 cout << "Segment: " << seg.edgenr << endl << "\t is REVERSE of the previous" << endl;
901 } // loop on geomEdge ancestors
903 if ( quadHelper ) // remember medium nodes of sub-meshes
905 SMDS_ElemIteratorPtr edges = smDS->GetElements();
906 while ( edges->more() )
908 const SMDS_MeshElement* e = edges->next();
909 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshEdge*>( e )))
915 } // case TopAbs_EDGE
917 case TopAbs_FACE: { // FACE
918 // ----------------------
919 const TopoDS_Face& geomFace = TopoDS::Face( sm->GetSubShape() );
920 helper.SetSubShape( geomFace );
921 bool isInternalFace = ( geomFace.Orientation() == TopAbs_INTERNAL );
923 // Find solids the geomFace bounds
924 int solidID1 = 0, solidID2 = 0;
925 StdMeshers_QuadToTriaAdaptor* quadAdaptor =
926 dynamic_cast<StdMeshers_QuadToTriaAdaptor*>( proxyMesh.get() );
929 solidID1 = occgeom.somap.FindIndex( quadAdaptor->GetShape() );
933 PShapeIteratorPtr solidIt = helper.GetAncestors( geomFace, *sm->GetFather(), TopAbs_SOLID);
934 while ( const TopoDS_Shape * solid = solidIt->next() )
936 int id = occgeom.somap.FindIndex ( *solid );
937 if ( solidID1 && id != solidID1 ) solidID2 = id;
941 // Add ng face descriptors of meshed faces
943 ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(faceNgID, solidID1, solidID2, 0));
945 // if second oreder is required, even already meshed faces must be passed to NETGEN
946 int fID = occgeom.fmap.Add( geomFace );
947 while ( fID < faceNgID ) // geomFace is already in occgeom.fmap, add a copy
948 fID = occgeom.fmap.Add( BRepBuilderAPI_Copy( geomFace, /*copyGeom=*/false ));
949 // Problem with the second order in a quadrangular mesh remains.
950 // 1) All quadrangles generated by NETGEN are moved to an inexistent face
951 // by FillSMesh() (find "AddFaceDescriptor")
952 // 2) Temporary triangles generated by StdMeshers_QuadToTriaAdaptor
953 // are on faces where quadrangles were.
954 // Due to these 2 points, wrong geom faces are used while conversion to qudratic
955 // of the mentioned above quadrangles and triangles
957 // Orient the face correctly in solidID1 (issue 0020206)
958 bool reverse = false;
960 TopoDS_Shape solid = occgeom.somap( solidID1 );
961 TopAbs_Orientation faceOriInSolid = helper.GetSubShapeOri( solid, geomFace );
962 if ( faceOriInSolid >= 0 )
964 helper.IsReversedSubMesh( TopoDS::Face( geomFace.Oriented( faceOriInSolid )));
967 // Add surface elements
969 netgen::Element2d tri(3);
970 tri.SetIndex ( faceNgID );
973 #ifdef DUMP_TRIANGLES
974 cout << "SMESH face " << helper.GetMeshDS()->ShapeToIndex( geomFace )
975 << " internal="<<isInternalFace << endl;
978 smDS = proxyMesh->GetSubMesh( geomFace );
980 SMDS_ElemIteratorPtr faces = smDS->GetElements();
981 while ( faces->more() )
983 const SMDS_MeshElement* f = faces->next();
984 if ( f->NbNodes() % 3 != 0 ) // not triangle
986 PShapeIteratorPtr solidIt=helper.GetAncestors(geomFace,*sm->GetFather(),TopAbs_SOLID);
987 if ( const TopoDS_Shape * solid = solidIt->next() )
988 sm = _mesh->GetSubMesh( *solid );
989 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
990 smError.reset( new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH,"Not triangle submesh"));
991 smError->myBadElements.push_back( f );
995 for ( int i = 0; i < 3; ++i )
997 const SMDS_MeshNode* node = f->GetNode( i ), * inFaceNode=0;
999 // get node UV on face
1000 int shapeID = node->getshapeId();
1001 if ( helper.IsSeamShape( shapeID ))
1002 if ( helper.IsSeamShape( f->GetNodeWrap( i+1 )->getshapeId() ))
1003 inFaceNode = f->GetNodeWrap( i-1 );
1005 inFaceNode = f->GetNodeWrap( i+1 );
1006 gp_XY uv = helper.GetNodeUV( geomFace, node, inFaceNode );
1008 int ind = reverse ? 3-i : i+1;
1009 tri.GeomInfoPi(ind).u = uv.X();
1010 tri.GeomInfoPi(ind).v = uv.Y();
1011 tri.PNum (ind) = ngNodeId( node, ngMesh, nodeNgIdMap );
1014 ngMesh.AddSurfaceElement (tri);
1015 #ifdef DUMP_TRIANGLES
1016 cout << tri << endl;
1019 if ( isInternalFace )
1021 swap( tri[1], tri[2] );
1022 ngMesh.AddSurfaceElement (tri);
1023 #ifdef DUMP_TRIANGLES
1024 cout << tri << endl;
1029 if ( quadHelper ) // remember medium nodes of sub-meshes
1031 SMDS_ElemIteratorPtr faces = smDS->GetElements();
1032 while ( faces->more() )
1034 const SMDS_MeshElement* f = faces->next();
1035 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshFace*>( f )))
1041 } // case TopAbs_FACE
1043 case TopAbs_VERTEX: { // VERTEX
1044 // --------------------------
1045 // issue 0021405. Add node only if a VERTEX is shared by a not meshed EDGE,
1046 // else netgen removes a free node and nodeVector becomes invalid
1047 PShapeIteratorPtr ansIt = helper.GetAncestors( sm->GetSubShape(),
1051 while ( const TopoDS_Shape* e = ansIt->next() )
1053 SMESH_subMesh* eSub = helper.GetMesh()->GetSubMesh( *e );
1054 if (( toAdd = ( eSub->IsEmpty() && !SMESH_Algo::isDegenerated( TopoDS::Edge( *e )))))
1059 SMDS_NodeIteratorPtr nodeIt = smDS->GetNodes();
1060 if ( nodeIt->more() )
1061 ngNodeId( nodeIt->next(), ngMesh, nodeNgIdMap );
1067 } // loop on submeshes
1070 nodeVec.resize( ngMesh.GetNP() + 1 );
1071 TNode2IdMap::iterator node_NgId, nodeNgIdEnd = nodeNgIdMap.end();
1072 for ( node_NgId = nodeNgIdMap.begin(); node_NgId != nodeNgIdEnd; ++node_NgId)
1073 nodeVec[ node_NgId->second ] = node_NgId->first;
1078 //================================================================================
1080 * \brief Duplicate mesh faces on internal geom faces
1082 //================================================================================
1084 void NETGENPlugin_Mesher::FixIntFaces(const netgen::OCCGeometry& occgeom,
1085 netgen::Mesh& ngMesh,
1086 NETGENPlugin_Internals& internalShapes)
1088 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1090 // find ng indices of internal faces
1092 for ( int ngFaceID = 1; ngFaceID <= occgeom.fmap.Extent(); ++ngFaceID )
1094 int smeshID = meshDS->ShapeToIndex( occgeom.fmap( ngFaceID ));
1095 if ( internalShapes.isInternalShape( smeshID ))
1096 ngFaceIds.insert( ngFaceID );
1098 if ( !ngFaceIds.empty() )
1101 int i, nbFaces = ngMesh.GetNSE();
1102 for (int i = 1; i <= nbFaces; ++i)
1104 netgen::Element2d elem = ngMesh.SurfaceElement(i);
1105 if ( ngFaceIds.count( elem.GetIndex() ))
1107 swap( elem[1], elem[2] );
1108 ngMesh.AddSurfaceElement (elem);
1114 //================================================================================
1116 * \brief Tries to heal the mesh on a FACE. The FACE is supposed to be partially
1117 * meshed due to NETGEN failure
1118 * \param [in] occgeom - geometry
1119 * \param [in,out] ngMesh - the mesh to fix
1120 * \param [inout] faceID - ID of the FACE to fix the mesh on
1121 * \return bool - is mesh is or becomes OK
1123 //================================================================================
1125 bool NETGENPlugin_Mesher::FixFaceMesh(const netgen::OCCGeometry& occgeom,
1126 netgen::Mesh& ngMesh,
1129 // we address a case where the FACE is almost fully meshed except small holes
1130 // of usually triangular shape at FACE boundary (IPAL52861)
1132 // The case appeared to be not simple: holes only look triangular but
1133 // indeed are a self intersecting polygon. A reason of the bug was in coincident
1134 // NG points on a seam edge. But the code below is very nice, leave it for
1139 if ( occgeom.fmap.Extent() < faceID )
1141 const TopoDS_Face& face = TopoDS::Face( occgeom.fmap( faceID ));
1143 // find free links on the FACE
1144 NCollection_Map<Link> linkMap;
1145 for ( int iF = 1; iF <= ngMesh.GetNSE(); ++iF )
1147 const netgen::Element2d& elem = ngMesh.SurfaceElement(iF);
1148 if ( faceID != elem.GetIndex() )
1150 int n0 = elem[ elem.GetNP() - 1 ];
1151 for ( int i = 0; i < elem.GetNP(); ++i )
1154 Link link( n0, n1 );
1155 if ( !linkMap.Add( link ))
1156 linkMap.Remove( link );
1160 // add/remove boundary links
1161 for ( int iSeg = 1; iSeg <= ngMesh.GetNSeg(); ++iSeg )
1163 const netgen::Segment& seg = ngMesh.LineSegment( iSeg );
1164 if ( seg.si != faceID ) // !edgeIDs.Contains( seg.edgenr ))
1166 Link link( seg[1], seg[0] ); // reverse!!!
1167 if ( !linkMap.Add( link ))
1168 linkMap.Remove( link );
1170 if ( linkMap.IsEmpty() )
1172 if ( linkMap.Extent() < 3 )
1175 // make triangles of the links
1177 netgen::Element2d tri(3);
1178 tri.SetIndex ( faceID );
1180 NCollection_Map<Link>::Iterator linkIt( linkMap );
1181 Link link1 = linkIt.Value();
1182 // look for a link connected to link1
1183 NCollection_Map<Link>::Iterator linkIt2 = linkIt;
1184 for ( linkIt2.Next(); linkIt2.More(); linkIt2.Next() )
1186 const Link& link2 = linkIt2.Value();
1187 if ( link2.IsConnected( link1 ))
1189 // look for a link connected to both link1 and link2
1190 NCollection_Map<Link>::Iterator linkIt3 = linkIt2;
1191 for ( linkIt3.Next(); linkIt3.More(); linkIt3.Next() )
1193 const Link& link3 = linkIt3.Value();
1194 if ( link3.IsConnected( link1 ) &&
1195 link3.IsConnected( link2 ) )
1200 tri[2] = ( link2.Contains( link1.n1 ) ? link2.n1 : link3.n1 );
1201 if ( tri[0] == tri[2] || tri[1] == tri[2] )
1203 ngMesh.AddSurfaceElement( tri );
1205 // prepare for the next tria search
1206 if ( linkMap.Extent() == 3 )
1208 linkMap.Remove( link3 );
1209 linkMap.Remove( link2 );
1211 linkMap.Remove( link1 );
1212 link1 = linkIt.Value();
1225 //================================================================================
1226 // define gp_XY_Subtracted pointer to function calling gp_XY::Subtracted(gp_XY)
1227 gp_XY_FunPtr(Subtracted);
1228 //gp_XY_FunPtr(Added);
1230 //================================================================================
1232 * \brief Evaluate distance between two 2d points along the surface
1234 //================================================================================
1236 double evalDist( const gp_XY& uv1,
1238 const Handle(Geom_Surface)& surf,
1239 const int stopHandler=-1)
1241 if ( stopHandler > 0 ) // continue recursion
1243 gp_XY mid = SMESH_MesherHelper::GetMiddleUV( surf, uv1, uv2 );
1244 return evalDist( uv1,mid, surf, stopHandler-1 ) + evalDist( mid,uv2, surf, stopHandler-1 );
1246 double dist3D = surf->Value( uv1.X(), uv1.Y() ).Distance( surf->Value( uv2.X(), uv2.Y() ));
1247 if ( stopHandler == 0 ) // stop recursion
1250 // start recursion if necessary
1251 double dist2D = SMESH_MesherHelper::ApplyIn2D(surf, uv1, uv2, gp_XY_Subtracted, 0).Modulus();
1252 if ( fabs( dist3D - dist2D ) < dist2D * 1e-10 )
1253 return dist3D; // equal parametrization of a planar surface
1255 return evalDist( uv1, uv2, surf, 3 ); // start recursion
1258 //================================================================================
1260 * \brief Data of vertex internal in geom face
1262 //================================================================================
1266 gp_XY uv; //!< UV in face parametric space
1267 int ngId; //!< ng id of corrsponding node
1268 gp_XY uvClose; //!< UV of closest boundary node
1269 int ngIdClose; //!< ng id of closest boundary node
1272 //================================================================================
1274 * \brief Data of vertex internal in solid
1276 //================================================================================
1280 int ngId; //!< ng id of corresponding node
1281 int ngIdClose; //!< ng id of closest 2d mesh element
1282 int ngIdCloseN; //!< ng id of closest node of the closest 2d mesh element
1285 inline double dist2(const netgen::MeshPoint& p1, const netgen::MeshPoint& p2)
1287 return gp_Pnt( NGPOINT_COORDS(p1)).SquareDistance( gp_Pnt( NGPOINT_COORDS(p2)));
1291 //================================================================================
1293 * \brief Make netgen take internal vertices in faces into account by adding
1294 * segments including internal vertices
1296 * This function works in supposition that 1D mesh is already computed in ngMesh
1298 //================================================================================
1300 void NETGENPlugin_Mesher::AddIntVerticesInFaces(const netgen::OCCGeometry& occgeom,
1301 netgen::Mesh& ngMesh,
1302 vector<const SMDS_MeshNode*>& nodeVec,
1303 NETGENPlugin_Internals& internalShapes)
1305 if ( nodeVec.size() < ngMesh.GetNP() )
1306 nodeVec.resize( ngMesh.GetNP(), 0 );
1308 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1309 SMESH_MesherHelper helper( internalShapes.getMesh() );
1311 const map<int,list<int> >& face2Vert = internalShapes.getFacesWithVertices();
1312 map<int,list<int> >::const_iterator f2v = face2Vert.begin();
1313 for ( ; f2v != face2Vert.end(); ++f2v )
1315 const TopoDS_Face& face = TopoDS::Face( meshDS->IndexToShape( f2v->first ));
1316 if ( face.IsNull() ) continue;
1317 int faceNgID = occgeom.fmap.FindIndex (face);
1318 if ( faceNgID < 0 ) continue;
1320 TopLoc_Location loc;
1321 Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc);
1323 helper.SetSubShape( face );
1324 helper.SetElementsOnShape( true );
1326 // Get data of internal vertices and add them to ngMesh
1328 multimap< double, TIntVData > dist2VData; // sort vertices by distance from boundary nodes
1330 int i, nbSegInit = ngMesh.GetNSeg();
1332 // boundary characteristics
1333 double totSegLen2D = 0;
1336 const list<int>& iVertices = f2v->second;
1337 list<int>::const_iterator iv = iVertices.begin();
1338 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1341 // get node on vertex
1342 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1343 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1346 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1347 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1348 nV = SMESH_Algo::VertexNode( V, meshDS );
1349 if ( !nV ) continue;
1352 netgen::MeshPoint mp( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1353 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1354 vData.ngId = ngMesh.GetNP();
1355 nodeVec.push_back( nV );
1359 vData.uv = helper.GetNodeUV( face, nV, 0, &uvOK );
1360 if ( !uvOK ) helper.CheckNodeUV( face, nV, vData.uv, BRep_Tool::Tolerance(V),/*force=*/1);
1362 // loop on all segments of the face to find the node closest to vertex and to count
1363 // average segment 2d length
1364 double closeDist2 = numeric_limits<double>::max(), dist2;
1366 for (i = 1; i <= ngMesh.GetNSeg(); ++i)
1368 netgen::Segment & seg = ngMesh.LineSegment(i);
1369 if ( seg.si != faceNgID ) continue;
1371 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1373 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1374 if ( ngIdLast == seg[ iEnd ] ) continue;
1375 dist2 = helper.ApplyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1376 if ( dist2 < closeDist2 )
1377 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1378 ngIdLast = seg[ iEnd ];
1382 totSegLen2D += helper.ApplyIn2D(surf, uv[0], uv[1], gp_XY_Subtracted, false).Modulus();
1386 dist2VData.insert( make_pair( closeDist2, vData ));
1389 if ( totNbSeg == 0 ) break;
1390 double avgSegLen2d = totSegLen2D / totNbSeg;
1392 // Loop on vertices to add segments
1394 multimap< double, TIntVData >::iterator dist_vData = dist2VData.begin();
1395 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1397 double closeDist2 = dist_vData->first, dist2;
1398 TIntVData & vData = dist_vData->second;
1400 // try to find more close node among segments added for internal vertices
1401 for (i = nbSegInit+1; i <= ngMesh.GetNSeg(); ++i)
1403 netgen::Segment & seg = ngMesh.LineSegment(i);
1404 if ( seg.si != faceNgID ) continue;
1406 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1408 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1409 dist2 = helper.ApplyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1410 if ( dist2 < closeDist2 )
1411 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1414 // decide whether to use the closest node as the second end of segment or to
1415 // create a new point
1416 int segEnd1 = vData.ngId;
1417 int segEnd2 = vData.ngIdClose; // to use closest node
1418 gp_XY uvV = vData.uv, uvP = vData.uvClose;
1419 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1420 double nodeDist2D = sqrt( closeDist2 );
1421 double nodeDist3D = evalDist( vData.uv, vData.uvClose, surf );
1422 bool avgLenOK = ( avgSegLen2d < 0.75 * nodeDist2D );
1423 bool hintLenOK = ( segLenHint < 0.75 * nodeDist3D );
1424 //cout << "uvV " << uvV.X() <<","<<uvV.Y() << " ";
1425 if ( hintLenOK || avgLenOK )
1427 // create a point between the closest node and V
1430 double r = min( 0.5, ( hintLenOK ? segLenHint/nodeDist3D : avgSegLen2d/nodeDist2D ));
1431 // direction from V to closet node in 2D
1432 gp_Dir2d v2n( helper.ApplyIn2D(surf, uvP, uvV, gp_XY_Subtracted, false ));
1434 uvP = vData.uv + r * nodeDist2D * v2n.XY();
1435 gp_Pnt P = surf->Value( uvP.X(), uvP.Y() ).Transformed( loc );
1437 netgen::MeshPoint mp( netgen::Point<3> (P.X(), P.Y(), P.Z()));
1438 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1439 segEnd2 = ngMesh.GetNP();
1440 //cout << "Middle " << r << " uv " << uvP.X() << "," << uvP.Y() << "( " << ngMesh.Point(segEnd2).X()<<","<<ngMesh.Point(segEnd2).Y()<<","<<ngMesh.Point(segEnd2).Z()<<" )"<< endl;
1441 SMDS_MeshNode * nP = helper.AddNode(P.X(), P.Y(), P.Z());
1442 nodeVec.push_back( nP );
1444 //else cout << "at Node " << " uv " << uvP.X() << "," << uvP.Y() << endl;
1447 netgen::Segment seg;
1449 if ( segEnd1 > segEnd2 ) swap( segEnd1, segEnd2 ), swap( uvV, uvP );
1450 seg[0] = segEnd1; // ng node id
1451 seg[1] = segEnd2; // ng node id
1452 seg.edgenr = ngMesh.GetNSeg() + 1;// segment id
1455 seg.epgeominfo[ 0 ].dist = 0; // param on curve
1456 seg.epgeominfo[ 0 ].u = uvV.X();
1457 seg.epgeominfo[ 0 ].v = uvV.Y();
1458 seg.epgeominfo[ 1 ].dist = 1; // param on curve
1459 seg.epgeominfo[ 1 ].u = uvP.X();
1460 seg.epgeominfo[ 1 ].v = uvP.Y();
1462 // seg.epgeominfo[ 0 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1463 // seg.epgeominfo[ 1 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1465 ngMesh.AddSegment (seg);
1467 // add reverse segment
1468 swap (seg[0], seg[1]);
1469 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1470 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1471 ngMesh.AddSegment (seg);
1477 //================================================================================
1479 * \brief Make netgen take internal vertices in solids into account by adding
1480 * faces including internal vertices
1482 * This function works in supposition that 2D mesh is already computed in ngMesh
1484 //================================================================================
1486 void NETGENPlugin_Mesher::AddIntVerticesInSolids(const netgen::OCCGeometry& occgeom,
1487 netgen::Mesh& ngMesh,
1488 vector<const SMDS_MeshNode*>& nodeVec,
1489 NETGENPlugin_Internals& internalShapes)
1491 #ifdef DUMP_TRIANGLES_SCRIPT
1492 // create a python script making a mesh containing triangles added for internal vertices
1493 ofstream py(DUMP_TRIANGLES_SCRIPT);
1494 py << "import SMESH"<< endl
1495 << "from salome.smesh import smeshBuilder"<<endl
1496 << "smesh = smeshBuilder.New(salome.myStudy)"<<endl
1497 << "m = smesh.Mesh(name='triangles')" << endl;
1499 if ( nodeVec.size() < ngMesh.GetNP() )
1500 nodeVec.resize( ngMesh.GetNP(), 0 );
1502 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1503 SMESH_MesherHelper helper( internalShapes.getMesh() );
1505 const map<int,list<int> >& so2Vert = internalShapes.getSolidsWithVertices();
1506 map<int,list<int> >::const_iterator s2v = so2Vert.begin();
1507 for ( ; s2v != so2Vert.end(); ++s2v )
1509 const TopoDS_Shape& solid = meshDS->IndexToShape( s2v->first );
1510 if ( solid.IsNull() ) continue;
1511 int solidNgID = occgeom.somap.FindIndex (solid);
1512 if ( solidNgID < 0 && !occgeom.somap.IsEmpty() ) continue;
1514 helper.SetSubShape( solid );
1515 helper.SetElementsOnShape( true );
1517 // find ng indices of faces within the solid
1519 for (TopExp_Explorer fExp(solid, TopAbs_FACE); fExp.More(); fExp.Next() )
1520 ngFaceIds.insert( occgeom.fmap.FindIndex( fExp.Current() ));
1521 if ( ngFaceIds.size() == 1 && *ngFaceIds.begin() == 0 )
1522 ngFaceIds.insert( 1 );
1524 // Get data of internal vertices and add them to ngMesh
1526 multimap< double, TIntVSoData > dist2VData; // sort vertices by distance from ng faces
1528 int i, nbFaceInit = ngMesh.GetNSE();
1530 // boundary characteristics
1531 double totSegLen = 0;
1534 const list<int>& iVertices = s2v->second;
1535 list<int>::const_iterator iv = iVertices.begin();
1536 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1539 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1541 // get node on vertex
1542 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1545 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1546 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1547 nV = SMESH_Algo::VertexNode( V, meshDS );
1548 if ( !nV ) continue;
1551 netgen::MeshPoint mpV( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1552 ngMesh.AddPoint ( mpV, 1, netgen::FIXEDPOINT );
1553 vData.ngId = ngMesh.GetNP();
1554 nodeVec.push_back( nV );
1556 // loop on all 2d elements to find the one closest to vertex and to count
1557 // average segment length
1558 double closeDist2 = numeric_limits<double>::max(), avgDist2;
1559 for (i = 1; i <= ngMesh.GetNSE(); ++i)
1561 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1562 if ( !ngFaceIds.count( elem.GetIndex() )) continue;
1564 multimap< double, int> dist2nID; // sort nodes of element by distance from V
1565 for ( int j = 0; j < elem.GetNP(); ++j)
1567 netgen::MeshPoint mp = ngMesh.Point( elem[j] );
1568 double d2 = dist2( mpV, mp );
1569 dist2nID.insert( make_pair( d2, elem[j] ));
1570 avgDist2 += d2 / elem.GetNP();
1572 totNbSeg++, totSegLen+= sqrt( dist2( mp, ngMesh.Point( elem[(j+1)%elem.GetNP()])));
1574 double dist = dist2nID.begin()->first; //avgDist2;
1575 if ( dist < closeDist2 )
1576 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= dist;
1578 dist2VData.insert( make_pair( closeDist2, vData ));
1581 if ( totNbSeg == 0 ) break;
1582 double avgSegLen = totSegLen / totNbSeg;
1584 // Loop on vertices to add triangles
1586 multimap< double, TIntVSoData >::iterator dist_vData = dist2VData.begin();
1587 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1589 double closeDist2 = dist_vData->first;
1590 TIntVSoData & vData = dist_vData->second;
1592 const netgen::MeshPoint& mpV = ngMesh.Point( vData.ngId );
1594 // try to find more close face among ones added for internal vertices
1595 for (i = nbFaceInit+1; i <= ngMesh.GetNSE(); ++i)
1597 double avgDist2 = 0;
1598 multimap< double, int> dist2nID;
1599 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1600 for ( int j = 0; j < elem.GetNP(); ++j)
1602 double d = dist2( mpV, ngMesh.Point( elem[j] ));
1603 dist2nID.insert( make_pair( d, elem[j] ));
1604 avgDist2 += d / elem.GetNP();
1605 if ( avgDist2 < closeDist2 )
1606 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= avgDist2;
1609 // sort nodes of the closest face by angle with vector from V to the closest node
1610 const double tol = numeric_limits<double>::min();
1611 map< double, int > angle2ID;
1612 const netgen::Element2d& closeFace = ngMesh.SurfaceElement( vData.ngIdClose );
1613 netgen::MeshPoint mp[2];
1614 mp[0] = ngMesh.Point( vData.ngIdCloseN );
1615 gp_XYZ p1( NGPOINT_COORDS( mp[0] ));
1616 gp_XYZ pV( NGPOINT_COORDS( mpV ));
1617 gp_Vec v2p1( pV, p1 );
1618 double distN1 = v2p1.Magnitude();
1619 if ( distN1 <= tol ) continue;
1621 for ( int j = 0; j < closeFace.GetNP(); ++j)
1623 mp[1] = ngMesh.Point( closeFace[j] );
1624 gp_Vec v2p( pV, gp_Pnt( NGPOINT_COORDS( mp[1] )) );
1625 angle2ID.insert( make_pair( v2p1.Angle( v2p ), closeFace[j]));
1627 // get node with angle of 60 degrees or greater
1628 map< double, int >::iterator angle_id = angle2ID.lower_bound( 60. * M_PI / 180. );
1629 if ( angle_id == angle2ID.end() ) angle_id = --angle2ID.end();
1630 const double minAngle = 30. * M_PI / 180.;
1631 const double angle = angle_id->first;
1632 bool angleOK = ( angle > minAngle );
1634 // find points to create a triangle
1635 netgen::Element2d tri(3);
1637 tri[0] = vData.ngId;
1638 tri[1] = vData.ngIdCloseN; // to use the closest nodes
1639 tri[2] = angle_id->second; // to use the node with best angle
1641 // decide whether to use the closest node and the node with best angle or to create new ones
1642 for ( int isBestAngleN = 0; isBestAngleN < 2; ++isBestAngleN )
1644 bool createNew = !angleOK, distOK = true;
1646 int triInd = isBestAngleN ? 2 : 1;
1647 mp[isBestAngleN] = ngMesh.Point( tri[triInd] );
1652 double distN2 = sqrt( dist2( mpV, mp[isBestAngleN]));
1653 createNew = ( fabs( distN2 - distN1 ) > 0.25 * distN1 );
1655 else if ( angle < tol )
1657 v2p1.SetX( v2p1.X() + 1e-3 );
1663 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1664 bool avgLenOK = ( avgSegLen < 0.75 * distN1 );
1665 bool hintLenOK = ( segLenHint < 0.75 * distN1 );
1666 createNew = (createNew || avgLenOK || hintLenOK );
1667 // we create a new node not closer than 0.5 to the closest face
1668 // in order not to clash with other close face
1669 double r = min( 0.5, ( hintLenOK ? segLenHint : avgSegLen ) / distN1 );
1670 distFromV = r * distN1;
1674 // create a new point, between the node and the vertex if angleOK
1675 gp_XYZ p( NGPOINT_COORDS( mp[isBestAngleN] ));
1676 gp_Vec v2p( pV, p ); v2p.Normalize();
1677 if ( isBestAngleN && !angleOK )
1678 p = p1 + gp_Dir( v2p.XYZ() - v2p1.XYZ()).XYZ() * distN1 * 0.95;
1680 p = pV + v2p.XYZ() * distFromV;
1682 if ( !isBestAngleN ) p1 = p, distN1 = distFromV;
1684 mp[isBestAngleN].SetPoint( netgen::Point<3> (p.X(), p.Y(), p.Z()));
1685 ngMesh.AddPoint ( mp[isBestAngleN], 1, netgen::SURFACEPOINT );
1686 tri[triInd] = ngMesh.GetNP();
1687 nodeVec.push_back( helper.AddNode( p.X(), p.Y(), p.Z()) );
1690 ngMesh.AddSurfaceElement (tri);
1691 swap( tri[1], tri[2] );
1692 ngMesh.AddSurfaceElement (tri);
1694 #ifdef DUMP_TRIANGLES_SCRIPT
1695 py << "n1 = m.AddNode( "<< mpV(0)<<", "<< mpV(1)<<", "<< mpV(2)<<") "<< endl
1696 << "n2 = m.AddNode( "<< mp[0](0)<<", "<< mp[0](1)<<", "<< mp[0](2)<<") "<< endl
1697 << "n3 = m.AddNode( "<< mp[1](0)<<", "<< mp[1](1)<<", "<< mp[1](2)<<" )" << endl
1698 << "m.AddFace([n1,n2,n3])" << endl;
1700 } // loop on internal vertices of a solid
1702 } // loop on solids with internal vertices
1705 //================================================================================
1707 * \brief Fill netgen mesh with segments of a FACE
1708 * \param ngMesh - netgen mesh
1709 * \param geom - container of OCCT geometry to mesh
1710 * \param wires - data of nodes on FACE boundary
1711 * \param helper - mesher helper holding the FACE
1712 * \param nodeVec - vector of nodes in which node index == netgen ID
1713 * \retval SMESH_ComputeErrorPtr - error description
1715 //================================================================================
1717 SMESH_ComputeErrorPtr
1718 NETGENPlugin_Mesher::AddSegmentsToMesh(netgen::Mesh& ngMesh,
1719 netgen::OCCGeometry& geom,
1720 const TSideVector& wires,
1721 SMESH_MesherHelper& helper,
1722 vector< const SMDS_MeshNode* > & nodeVec,
1723 const bool overrideMinH)
1725 // ----------------------------
1726 // Check wires and count nodes
1727 // ----------------------------
1729 for ( int iW = 0; iW < wires.size(); ++iW )
1731 StdMeshers_FaceSidePtr wire = wires[ iW ];
1732 if ( wire->MissVertexNode() )
1734 // Commented for issue 0020960. It worked for the case, let's wait for case where it doesn't.
1735 // It seems that there is no reason for this limitation
1737 // (new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH, "Missing nodes on vertices"));
1739 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1740 if ( uvPtVec.size() != wire->NbPoints() )
1741 return SMESH_ComputeError::New(COMPERR_BAD_INPUT_MESH,
1742 SMESH_Comment("Unexpected nb of points on wire ") << iW
1743 << ": " << uvPtVec.size()<<" != "<<wire->NbPoints());
1744 nbNodes += wire->NbPoints();
1746 nodeVec.reserve( nodeVec.size() + nbNodes + 1 );
1747 if ( nodeVec.empty() )
1748 nodeVec.push_back( 0 );
1750 // -----------------
1752 // -----------------
1754 const bool wasNgMeshEmpty = ( ngMesh.GetNP() < 1 ); /* true => this method is called by
1755 NETGENPlugin_NETGEN_2D_ONLY */
1757 // map for nodes on vertices since they can be shared between wires
1758 // ( issue 0020676, face_int_box.brep) and nodes built by NETGEN
1759 map<const SMDS_MeshNode*, int > node2ngID;
1760 if ( !wasNgMeshEmpty ) // fill node2ngID with nodes built by NETGEN
1762 set< int > subIDs; // ids of sub-shapes of the FACE
1763 for ( int iW = 0; iW < wires.size(); ++iW )
1765 StdMeshers_FaceSidePtr wire = wires[ iW ];
1766 for ( int iE = 0, nbE = wire->NbEdges(); iE < nbE; ++iE )
1768 subIDs.insert( wire->EdgeID( iE ));
1769 subIDs.insert( helper.GetMeshDS()->ShapeToIndex( wire->FirstVertex( iE )));
1772 for ( size_t ngID = 1; ngID < nodeVec.size(); ++ngID )
1773 if ( subIDs.count( nodeVec[ngID]->getshapeId() ))
1774 node2ngID.insert( make_pair( nodeVec[ngID], ngID ));
1777 const int solidID = 0, faceID = geom.fmap.FindIndex( helper.GetSubShape() );
1778 if ( ngMesh.GetNFD() < 1 )
1779 ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(faceID, solidID, solidID, 0));
1781 for ( int iW = 0; iW < wires.size(); ++iW )
1783 StdMeshers_FaceSidePtr wire = wires[ iW ];
1784 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1785 const int nbSegments = wire->NbPoints() - 1;
1787 // assure the 1st node to be in node2ngID, which is needed to correctly
1788 // "close chain of segments" (see below) in case if the 1st node is not
1789 // onVertex because it is on a Viscous layer
1790 node2ngID.insert( make_pair( uvPtVec[ 0 ].node, ngMesh.GetNP() + 1 ));
1792 // compute length of every segment
1793 vector<double> segLen( nbSegments );
1794 for ( int i = 0; i < nbSegments; ++i )
1795 segLen[i] = SMESH_TNodeXYZ( uvPtVec[ i ].node ).Distance( uvPtVec[ i+1 ].node );
1797 int edgeID = 1, posID = -2;
1798 bool isInternalWire = false;
1799 double vertexNormPar = 0;
1800 const int prevNbNGSeg = ngMesh.GetNSeg();
1801 for ( int i = 0; i < nbSegments; ++i ) // loop on segments
1803 // Add the first point of a segment
1805 const SMDS_MeshNode * n = uvPtVec[ i ].node;
1806 const int posShapeID = n->getshapeId();
1807 bool onVertex = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX );
1808 bool onEdge = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_EDGE );
1810 // skip nodes on degenerated edges
1811 if ( helper.IsDegenShape( posShapeID ) &&
1812 helper.IsDegenShape( uvPtVec[ i+1 ].node->getshapeId() ))
1815 int ngID1 = ngMesh.GetNP() + 1, ngID2 = ngID1+1;
1816 if ( onVertex || ( !wasNgMeshEmpty && onEdge ) || helper.IsRealSeam( posShapeID ))
1817 ngID1 = node2ngID.insert( make_pair( n, ngID1 )).first->second;
1818 if ( ngID1 > ngMesh.GetNP() )
1820 netgen::MeshPoint mp( netgen::Point<3> (n->X(), n->Y(), n->Z()) );
1821 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1822 nodeVec.push_back( n );
1824 else // n is in ngMesh already, and ngID2 in prev segment is wrong
1826 ngID2 = ngMesh.GetNP() + 1;
1827 if ( i > 0 ) // prev segment belongs to same wire
1829 netgen::Segment& prevSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
1836 netgen::Segment seg;
1838 seg[0] = ngID1; // ng node id
1839 seg[1] = ngID2; // ng node id
1840 seg.edgenr = ngMesh.GetNSeg() + 1; // ng segment id
1841 seg.si = faceID; // = geom.fmap.FindIndex (face);
1843 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1845 const UVPtStruct& pnt = uvPtVec[ i + iEnd ];
1847 seg.epgeominfo[ iEnd ].dist = pnt.param; // param on curve
1848 seg.epgeominfo[ iEnd ].u = pnt.u;
1849 seg.epgeominfo[ iEnd ].v = pnt.v;
1851 // find out edge id and node parameter on edge
1852 onVertex = ( pnt.normParam + 1e-10 > vertexNormPar );
1853 if ( onVertex || posShapeID != posID )
1856 double normParam = pnt.normParam;
1858 normParam = 0.5 * ( uvPtVec[ i ].normParam + uvPtVec[ i+1 ].normParam );
1859 int edgeIndexInWire = wire->EdgeIndex( normParam );
1860 vertexNormPar = wire->LastParameter( edgeIndexInWire );
1861 const TopoDS_Edge& edge = wire->Edge( edgeIndexInWire );
1862 edgeID = geom.emap.FindIndex( edge );
1864 isInternalWire = ( edge.Orientation() == TopAbs_INTERNAL );
1865 // if ( onVertex ) // param on curve is different on each of two edges
1866 // seg.epgeominfo[ iEnd ].dist = helper.GetNodeU( edge, pnt.node );
1868 seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
1871 ngMesh.AddSegment (seg);
1873 // restrict size of elements near the segment
1874 SMESH_TNodeXYZ np1( n ), np2( uvPtVec[ i+1 ].node );
1875 // get an average size of adjacent segments to avoid sharp change of
1876 // element size (regression on issue 0020452, note 0010898)
1877 int iPrev = SMESH_MesherHelper::WrapIndex( i-1, nbSegments );
1878 int iNext = SMESH_MesherHelper::WrapIndex( i+1, nbSegments );
1879 double sumH = segLen[ iPrev ] + segLen[ i ] + segLen[ iNext ];
1880 int nbSeg = ( int( segLen[ iPrev ] > sumH / 100.) +
1881 int( segLen[ i ] > sumH / 100.) +
1882 int( segLen[ iNext ] > sumH / 100.));
1884 RestrictLocalSize( ngMesh, 0.5*(np1+np2), sumH / nbSeg, overrideMinH );
1886 if ( isInternalWire )
1888 swap (seg[0], seg[1]);
1889 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1890 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1891 ngMesh.AddSegment (seg);
1893 } // loop on segments on a wire
1895 // close chain of segments
1896 if ( nbSegments > 0 )
1898 netgen::Segment& lastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() - int( isInternalWire));
1899 const SMDS_MeshNode * lastNode = uvPtVec.back().node;
1900 lastSeg[1] = node2ngID.insert( make_pair( lastNode, lastSeg[1] )).first->second;
1901 if ( lastSeg[1] > ngMesh.GetNP() )
1903 netgen::MeshPoint mp( netgen::Point<3> (lastNode->X(), lastNode->Y(), lastNode->Z()) );
1904 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1905 nodeVec.push_back( lastNode );
1907 if ( isInternalWire )
1909 netgen::Segment& realLastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
1910 realLastSeg[0] = lastSeg[1];
1914 #ifdef DUMP_SEGMENTS
1915 cout << "BEGIN WIRE " << iW << endl;
1916 for ( int i = prevNbNGSeg+1; i <= ngMesh.GetNSeg(); ++i )
1918 netgen::Segment& seg = ngMesh.LineSegment( i );
1920 netgen::Segment& prevSeg = ngMesh.LineSegment( i-1 );
1921 if ( seg[0] == prevSeg[1] && seg[1] == prevSeg[0] )
1923 cout << "Segment: " << seg.edgenr << endl << "\tis REVRESE of the previous one" << endl;
1927 cout << "Segment: " << seg.edgenr << endl
1928 << "\tp1: " << seg[0] << " n" << nodeVec[ seg[0]]->GetID() << endl
1929 << "\tp2: " << seg[1] << " n" << nodeVec[ seg[1]]->GetID() << endl
1930 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
1931 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
1932 << "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
1933 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
1934 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl
1935 << "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
1937 cout << "--END WIRE " << iW << endl;
1940 } // loop on WIREs of a FACE
1942 // add a segment instead of an internal vertex
1943 if ( wasNgMeshEmpty )
1945 NETGENPlugin_Internals intShapes( *helper.GetMesh(), helper.GetSubShape(), /*is3D=*/false );
1946 AddIntVerticesInFaces( geom, ngMesh, nodeVec, intShapes );
1948 ngMesh.CalcSurfacesOfNode();
1953 //================================================================================
1955 * \brief Fill SMESH mesh according to contents of netgen mesh
1956 * \param occgeo - container of OCCT geometry to mesh
1957 * \param ngMesh - netgen mesh
1958 * \param initState - bn of entities in netgen mesh before computing
1959 * \param sMesh - SMESH mesh to fill in
1960 * \param nodeVec - vector of nodes in which node index == netgen ID
1961 * \param comment - returns problem description
1962 * \param quadHelper - holder of medium nodes of sub-meshes
1963 * \retval int - error
1965 //================================================================================
1967 int NETGENPlugin_Mesher::FillSMesh(const netgen::OCCGeometry& occgeo,
1968 netgen::Mesh& ngMesh,
1969 const NETGENPlugin_ngMeshInfo& initState,
1971 std::vector<const SMDS_MeshNode*>& nodeVec,
1972 SMESH_Comment& comment,
1973 SMESH_MesherHelper* quadHelper)
1975 int nbNod = ngMesh.GetNP();
1976 int nbSeg = ngMesh.GetNSeg();
1977 int nbFac = ngMesh.GetNSE();
1978 int nbVol = ngMesh.GetNE();
1980 SMESHDS_Mesh* meshDS = sMesh.GetMeshDS();
1982 // quadHelper is used for either
1983 // 1) making quadratic elements when a lower dimention mesh is loaded
1984 // to SMESH before convertion to quadratic by NETGEN
1985 // 2) sewing of quadratic elements with quadratic elements of sub-meshes
1986 if ( quadHelper && !quadHelper->GetIsQuadratic() && quadHelper->GetTLinkNodeMap().empty() )
1989 // -------------------------------------
1990 // Create and insert nodes into nodeVec
1991 // -------------------------------------
1993 nodeVec.resize( nbNod + 1 );
1994 int i, nbInitNod = initState._nbNodes;
1995 for (i = nbInitNod+1; i <= nbNod; ++i )
1997 const netgen::MeshPoint& ngPoint = ngMesh.Point(i);
1998 SMDS_MeshNode* node = NULL;
1999 TopoDS_Vertex aVert;
2000 // First, netgen creates nodes on vertices in occgeo.vmap,
2001 // so node index corresponds to vertex index
2002 // but (issue 0020776) netgen does not create nodes with equal coordinates
2003 if ( i-nbInitNod <= occgeo.vmap.Extent() )
2005 gp_Pnt p ( NGPOINT_COORDS(ngPoint) );
2006 for (int iV = i-nbInitNod; aVert.IsNull() && iV <= occgeo.vmap.Extent(); ++iV)
2008 aVert = TopoDS::Vertex( occgeo.vmap( iV ));
2009 gp_Pnt pV = BRep_Tool::Pnt( aVert );
2010 if ( p.SquareDistance( pV ) > 1e-20 )
2013 node = const_cast<SMDS_MeshNode*>( SMESH_Algo::VertexNode( aVert, meshDS ));
2016 if (!node) // node not found on vertex
2018 node = meshDS->AddNode( NGPOINT_COORDS( ngPoint ));
2019 if (!aVert.IsNull())
2020 meshDS->SetNodeOnVertex(node, aVert);
2025 // -------------------------------------------
2026 // Create mesh segments along geometric edges
2027 // -------------------------------------------
2029 int nbInitSeg = initState._nbSegments;
2030 for (i = nbInitSeg+1; i <= nbSeg; ++i )
2032 const netgen::Segment& seg = ngMesh.LineSegment(i);
2034 int pinds[3] = { seg.pnums[0], seg.pnums[1], seg.pnums[2] };
2037 for (int j=0; j < 3; ++j)
2039 int pind = pinds[j];
2040 if (pind <= 0 || !nodeVec_ACCESS(pind))
2048 int aGeomEdgeInd = seg.epgeominfo[j].edgenr;
2049 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
2050 aEdge = TopoDS::Edge(occgeo.emap(aGeomEdgeInd));
2052 param = seg.epgeominfo[j].dist;
2055 else // middle point
2057 param = param2 * 0.5;
2059 if (!aEdge.IsNull() && nodeVec_ACCESS(pind)->getshapeId() < 1)
2061 meshDS->SetNodeOnEdge(nodeVec_ACCESS(pind), aEdge, param);
2066 SMDS_MeshEdge* edge = 0;
2067 if (nbp == 2) // second order ?
2069 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1])))
2071 if ( quadHelper ) // final mesh must be quadratic
2072 edge = quadHelper->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
2074 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
2078 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
2079 nodeVec_ACCESS(pinds[2])))
2081 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
2082 nodeVec_ACCESS(pinds[2]));
2086 if ( comment.empty() ) comment << "Cannot create a mesh edge";
2087 MESSAGE("Cannot create a mesh edge");
2088 nbSeg = nbFac = nbVol = 0;
2091 if ( !aEdge.IsNull() && edge->getshapeId() < 1 )
2092 meshDS->SetMeshElementOnShape(edge, aEdge);
2094 else if ( comment.empty() )
2096 comment << "Invalid netgen segment #" << i;
2100 // ----------------------------------------
2101 // Create mesh faces along geometric faces
2102 // ----------------------------------------
2104 int nbInitFac = initState._nbFaces;
2105 int quadFaceID = ngMesh.GetNFD() + 1;
2106 if ( nbInitFac < nbFac )
2107 // add a faces descriptor to exclude qudrangle elements generated by NETGEN
2108 // from computation of 3D mesh
2109 ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(quadFaceID, /*solid1=*/0, /*solid2=*/0, 0));
2111 vector<const SMDS_MeshNode*> nodes;
2112 for (i = nbInitFac+1; i <= nbFac; ++i )
2114 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
2115 int aGeomFaceInd = elem.GetIndex();
2117 if (aGeomFaceInd > 0 && aGeomFaceInd <= occgeo.fmap.Extent())
2118 aFace = TopoDS::Face(occgeo.fmap(aGeomFaceInd));
2120 for (int j=1; j <= elem.GetNP(); ++j)
2122 int pind = elem.PNum(j);
2123 if ( pind < 1 || pind >= nodeVec.size() )
2125 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind))
2127 nodes.push_back( node );
2128 if (!aFace.IsNull() && node->getshapeId() < 1)
2130 const netgen::PointGeomInfo& pgi = elem.GeomInfoPi(j);
2131 meshDS->SetNodeOnFace(node, aFace, pgi.u, pgi.v);
2135 if ( nodes.size() != elem.GetNP() )
2137 if ( comment.empty() )
2138 comment << "Invalid netgen 2d element #" << i;
2139 continue; // bad node ids
2141 SMDS_MeshFace* face = NULL;
2142 switch (elem.GetType())
2145 if ( quadHelper ) // final mesh must be quadratic
2146 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2]);
2148 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2]);
2151 if ( quadHelper ) // final mesh must be quadratic
2152 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2154 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2155 // exclude qudrangle elements from computation of 3D mesh
2156 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2159 nodes[5] = mediumNode( nodes[0],nodes[1],nodes[5], quadHelper );
2160 nodes[3] = mediumNode( nodes[1],nodes[2],nodes[3], quadHelper );
2161 nodes[4] = mediumNode( nodes[2],nodes[0],nodes[4], quadHelper );
2162 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[5],nodes[3],nodes[4]);
2165 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2166 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2167 nodes[5] = mediumNode( nodes[2],nodes[3],nodes[5], quadHelper );
2168 nodes[6] = mediumNode( nodes[3],nodes[0],nodes[6], quadHelper );
2169 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3],
2170 nodes[4],nodes[7],nodes[5],nodes[6]);
2171 // exclude qudrangle elements from computation of 3D mesh
2172 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2175 MESSAGE("NETGEN created a face of unexpected type, ignoring");
2180 if ( comment.empty() ) comment << "Cannot create a mesh face";
2181 MESSAGE("Cannot create a mesh face");
2182 nbSeg = nbFac = nbVol = 0;
2185 if (!aFace.IsNull())
2186 meshDS->SetMeshElementOnShape(face, aFace);
2189 // ------------------
2190 // Create tetrahedra
2191 // ------------------
2193 for (i = 1; i <= nbVol; ++i)
2195 const netgen::Element& elem = ngMesh.VolumeElement(i);
2196 int aSolidInd = elem.GetIndex();
2197 TopoDS_Solid aSolid;
2198 if (aSolidInd > 0 && aSolidInd <= occgeo.somap.Extent())
2199 aSolid = TopoDS::Solid(occgeo.somap(aSolidInd));
2201 for (int j=1; j <= elem.GetNP(); ++j)
2203 int pind = elem.PNum(j);
2204 if ( pind < 1 || pind >= nodeVec.size() )
2206 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind) )
2208 nodes.push_back(node);
2209 if ( !aSolid.IsNull() && node->getshapeId() < 1 )
2210 meshDS->SetNodeInVolume(node, aSolid);
2213 if ( nodes.size() != elem.GetNP() )
2215 if ( comment.empty() )
2216 comment << "Invalid netgen 3d element #" << i;
2219 SMDS_MeshVolume* vol = NULL;
2220 switch (elem.GetType())
2223 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3]);
2226 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2227 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2228 nodes[5] = mediumNode( nodes[2],nodes[0],nodes[5], quadHelper );
2229 nodes[6] = mediumNode( nodes[0],nodes[3],nodes[6], quadHelper );
2230 nodes[8] = mediumNode( nodes[1],nodes[3],nodes[8], quadHelper );
2231 nodes[9] = mediumNode( nodes[2],nodes[3],nodes[9], quadHelper );
2232 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3],
2233 nodes[4],nodes[7],nodes[5],nodes[6],nodes[8],nodes[9]);
2236 MESSAGE("NETGEN created a volume of unexpected type, ignoring");
2241 if ( comment.empty() ) comment << "Cannot create a mesh volume";
2242 MESSAGE("Cannot create a mesh volume");
2243 nbSeg = nbFac = nbVol = 0;
2246 if (!aSolid.IsNull())
2247 meshDS->SetMeshElementOnShape(vol, aSolid);
2249 return comment.empty() ? 0 : 1;
2254 //================================================================================
2256 * \brief Restrict size of elements on the given edge
2258 //================================================================================
2260 void setLocalSize(const TopoDS_Edge& edge,
2264 if ( size <= std::numeric_limits<double>::min() )
2266 Standard_Real u1, u2;
2267 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, u1, u2);
2268 if ( curve.IsNull() )
2270 TopoDS_Iterator vIt( edge );
2271 if ( !vIt.More() ) return;
2272 gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex( vIt.Value() ));
2273 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
2277 const int nb = (int)( 1.5 * SMESH_Algo::EdgeLength( edge ) / size );
2278 Standard_Real delta = (u2-u1)/nb;
2279 for(int i=0; i<nb; i++)
2281 Standard_Real u = u1 + delta*i;
2282 gp_Pnt p = curve->Value(u);
2283 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
2284 netgen::Point3d pi(p.X(), p.Y(), p.Z());
2285 double resultSize = mesh.GetH(pi);
2286 if ( resultSize - size > 0.1*size )
2287 // netgen does restriction iff oldH/newH > 1.2 (localh.cpp:136)
2288 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), resultSize/1.201 );
2293 //================================================================================
2295 * \brief Convert error into text
2297 //================================================================================
2299 std::string text(int err)
2304 SMESH_Comment("Error in netgen::OCCGenerateMesh() at ") << netgen::multithread.task;
2307 //================================================================================
2309 * \brief Convert exception into text
2311 //================================================================================
2313 std::string text(Standard_Failure& ex)
2315 SMESH_Comment str("Exception in netgen::OCCGenerateMesh()");
2316 str << " at " << netgen::multithread.task
2317 << ": " << ex.DynamicType()->Name();
2318 if ( ex.GetMessageString() && strlen( ex.GetMessageString() ))
2319 str << ": " << ex.GetMessageString();
2322 //================================================================================
2324 * \brief Convert exception into text
2326 //================================================================================
2328 std::string text(netgen::NgException& ex)
2330 SMESH_Comment str("NgException");
2331 if ( strlen( netgen::multithread.task ) > 0 )
2332 str << " at " << netgen::multithread.task;
2333 str << ": " << ex.What();
2337 //================================================================================
2339 * \brief Looks for triangles lying on a SOLID
2341 //================================================================================
2343 bool hasBadElemOnSolid( const list<const SMDS_MeshElement*>& elems,
2344 SMESH_subMesh* solidSM )
2346 TopTools_IndexedMapOfShape solidSubs;
2347 TopExp::MapShapes( solidSM->GetSubShape(), solidSubs );
2348 SMESHDS_Mesh* mesh = solidSM->GetFather()->GetMeshDS();
2350 list<const SMDS_MeshElement*>::const_iterator e = elems.begin();
2351 for ( ; e != elems.end(); ++e )
2353 const SMDS_MeshElement* elem = *e;
2354 // if ( elem->GetType() != SMDSAbs_Face ) -- 23047
2356 int nbNodesOnSolid = 0, nbNodes = elem->NbNodes();
2357 SMDS_NodeIteratorPtr nIt = elem->nodeIterator();
2358 while ( nIt->more() )
2360 const SMDS_MeshNode* n = nIt->next();
2361 const TopoDS_Shape& s = mesh->IndexToShape( n->getshapeId() );
2362 nbNodesOnSolid += ( !s.IsNull() && solidSubs.Contains( s ));
2363 if ( nbNodesOnSolid > 2 ||
2364 nbNodesOnSolid == nbNodes)
2371 const double edgeMeshingTime = 0.001;
2372 const double faceMeshingTime = 0.019;
2373 const double edgeFaceMeshingTime = edgeMeshingTime + faceMeshingTime;
2374 const double faceOptimizTime = 0.06;
2375 const double voluMeshingTime = 0.15;
2376 const double volOptimizeTime = 0.77;
2379 //=============================================================================
2381 * Here we are going to use the NETGEN mesher
2383 //=============================================================================
2385 bool NETGENPlugin_Mesher::Compute()
2387 NETGENPlugin_NetgenLibWrapper ngLib;
2389 netgen::MeshingParameters& mparams = netgen::mparam;
2390 MESSAGE("Compute with:\n"
2391 " max size = " << mparams.maxh << "\n"
2392 " segments per edge = " << mparams.segmentsperedge);
2394 " growth rate = " << mparams.grading << "\n"
2395 " elements per radius = " << mparams.curvaturesafety << "\n"
2396 " second order = " << mparams.secondorder << "\n"
2397 " quad allowed = " << mparams.quad << "\n"
2398 " surface curvature = " << mparams.uselocalh << "\n"
2399 " fuse edges = " << netgen::merge_solids);
2401 SMESH_ComputeErrorPtr error = SMESH_ComputeError::New();
2402 SMESH_MesherHelper quadHelper( *_mesh );
2403 quadHelper.SetIsQuadratic( mparams.secondorder );
2405 static string debugFile = "/tmp/ngMesh.py"; /* to call toPython( _ngMesh, debugFile )
2406 while debugging netgen */
2407 // -------------------------
2408 // Prepare OCC geometry
2409 // -------------------------
2411 netgen::OCCGeometry occgeo;
2412 list< SMESH_subMesh* > meshedSM[3]; // for 0-2 dimensions
2413 NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
2414 PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
2417 _totalTime = edgeFaceMeshingTime;
2419 _totalTime += faceOptimizTime;
2421 _totalTime += voluMeshingTime + ( _optimize ? volOptimizeTime : 0 );
2422 double doneTime = 0;
2425 _curShapeIndex = -1;
2427 // -------------------------
2428 // Generate the mesh
2429 // -------------------------
2432 NETGENPlugin_ngMeshInfo initState; // it remembers size of ng mesh equal to size of Smesh
2434 SMESH_Comment comment;
2437 // vector of nodes in which node index == netgen ID
2438 vector< const SMDS_MeshNode* > nodeVec;
2446 // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
2447 mparams.uselocalh = false;
2448 mparams.grading = 0.8; // not limitited size growth
2450 if ( _simpleHyp->GetNumberOfSegments() )
2452 mparams.maxh = occgeo.boundingbox.Diam();
2455 mparams.maxh = _simpleHyp->GetLocalLength();
2458 if ( mparams.maxh == 0.0 )
2459 mparams.maxh = occgeo.boundingbox.Diam();
2460 if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
2461 mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
2463 // Local size on faces
2464 occgeo.face_maxh = mparams.maxh;
2466 // Let netgen create _ngMesh and calculate element size on not meshed shapes
2470 int startWith = netgen::MESHCONST_ANALYSE;
2471 int endWith = netgen::MESHCONST_ANALYSE;
2476 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2478 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2480 if(netgen::multithread.terminate)
2483 comment << text(err);
2485 catch (Standard_Failure& ex)
2487 comment << text(ex);
2489 err = 0; //- MESHCONST_ANALYSE isn't so important step
2492 ngLib.setMesh(( Ng_Mesh*) _ngMesh );
2494 _ngMesh->ClearFaceDescriptors(); // we make descriptors our-self
2498 // Pass 1D simple parameters to NETGEN
2499 // --------------------------------
2500 int nbSeg = _simpleHyp->GetNumberOfSegments();
2501 double segSize = _simpleHyp->GetLocalLength();
2502 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
2504 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
2506 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
2507 setLocalSize( e, segSize, *_ngMesh );
2510 else // if ( ! _simpleHyp )
2512 // Local size on vertices and edges
2513 // --------------------------------
2514 for(std::map<int,double>::const_iterator it=EdgeId2LocalSize.begin(); it!=EdgeId2LocalSize.end(); it++)
2516 int key = (*it).first;
2517 double hi = (*it).second;
2518 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2519 const TopoDS_Edge& e = TopoDS::Edge(shape);
2520 setLocalSize( e, hi, *_ngMesh );
2522 for(std::map<int,double>::const_iterator it=VertexId2LocalSize.begin(); it!=VertexId2LocalSize.end(); it++)
2524 int key = (*it).first;
2525 double hi = (*it).second;
2526 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2527 const TopoDS_Vertex& v = TopoDS::Vertex(shape);
2528 gp_Pnt p = BRep_Tool::Pnt(v);
2529 NETGENPlugin_Mesher::RestrictLocalSize( *_ngMesh, p.XYZ(), hi );
2531 for(map<int,double>::const_iterator it=FaceId2LocalSize.begin();
2532 it!=FaceId2LocalSize.end(); it++)
2534 int key = (*it).first;
2535 double val = (*it).second;
2536 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2537 int faceNgID = occgeo.fmap.FindIndex(shape);
2538 occgeo.SetFaceMaxH(faceNgID, val);
2539 for ( TopExp_Explorer edgeExp( shape, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
2540 setLocalSize( TopoDS::Edge( edgeExp.Current() ), val, *_ngMesh );
2544 // Precompute internal edges (issue 0020676) in order to
2545 // add mesh on them correctly (twice) to netgen mesh
2546 if ( !err && internals.hasInternalEdges() )
2548 // load internal shapes into OCCGeometry
2549 netgen::OCCGeometry intOccgeo;
2550 internals.getInternalEdges( intOccgeo.fmap, intOccgeo.emap, intOccgeo.vmap, meshedSM );
2551 intOccgeo.boundingbox = occgeo.boundingbox;
2552 intOccgeo.shape = occgeo.shape;
2553 intOccgeo.face_maxh.SetSize(intOccgeo.fmap.Extent());
2554 intOccgeo.face_maxh = netgen::mparam.maxh;
2555 netgen::Mesh *tmpNgMesh = NULL;
2559 // compute local H on internal shapes in the main mesh
2560 //OCCSetLocalMeshSize(intOccgeo, *_ngMesh); it deletes _ngMesh->localH
2562 // let netgen create a temporary mesh
2564 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2566 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2568 if(netgen::multithread.terminate)
2571 // copy LocalH from the main to temporary mesh
2572 initState.transferLocalH( _ngMesh, tmpNgMesh );
2574 // compute mesh on internal edges
2575 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2577 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2579 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2581 comment << text(err);
2583 catch (Standard_Failure& ex)
2585 comment << text(ex);
2588 initState.restoreLocalH( tmpNgMesh );
2590 // fill SMESH by netgen mesh
2591 vector< const SMDS_MeshNode* > tmpNodeVec;
2592 FillSMesh( intOccgeo, *tmpNgMesh, initState, *_mesh, tmpNodeVec, comment );
2593 err = ( err || !comment.empty() );
2595 nglib::Ng_DeleteMesh((nglib::Ng_Mesh*)tmpNgMesh);
2598 // Fill _ngMesh with nodes and segments of computed submeshes
2601 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_0D ]) &&
2602 FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_1D ], &quadHelper));
2604 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2609 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2614 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2616 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2618 if(netgen::multithread.terminate)
2621 comment << text(err);
2623 catch (Standard_Failure& ex)
2625 comment << text(ex);
2630 _ticTime = ( doneTime += edgeMeshingTime ) / _totalTime / _progressTic;
2632 mparams.uselocalh = true; // restore as it is used at surface optimization
2634 // ---------------------
2635 // compute surface mesh
2636 // ---------------------
2639 // Pass 2D simple parameters to NETGEN
2641 if ( double area = _simpleHyp->GetMaxElementArea() ) {
2643 mparams.maxh = sqrt(2. * area/sqrt(3.0));
2644 mparams.grading = 0.4; // moderate size growth
2647 // length from edges
2648 if ( _ngMesh->GetNSeg() ) {
2649 double edgeLength = 0;
2650 TopTools_MapOfShape visitedEdges;
2651 for ( TopExp_Explorer exp( _shape, TopAbs_EDGE ); exp.More(); exp.Next() )
2652 if( visitedEdges.Add(exp.Current()) )
2653 edgeLength += SMESH_Algo::EdgeLength( TopoDS::Edge( exp.Current() ));
2654 // we have to multiply length by 2 since for each TopoDS_Edge there
2655 // are double set of NETGEN edges, in other words, we have to
2656 // divide _ngMesh->GetNSeg() by 2.
2657 mparams.maxh = 2*edgeLength / _ngMesh->GetNSeg();
2660 mparams.maxh = 1000;
2662 mparams.grading = 0.2; // slow size growth
2664 mparams.quad = _simpleHyp->GetAllowQuadrangles();
2665 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2666 _ngMesh->SetGlobalH (mparams.maxh);
2667 netgen::Box<3> bb = occgeo.GetBoundingBox();
2668 bb.Increase (bb.Diam()/20);
2669 _ngMesh->SetLocalH (bb.PMin(), bb.PMax(), mparams.grading);
2672 // Care of vertices internal in faces (issue 0020676)
2673 if ( internals.hasInternalVertexInFace() )
2675 // store computed segments in SMESH in order not to create SMESH
2676 // edges for ng segments added by AddIntVerticesInFaces()
2677 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2678 // add segments to faces with internal vertices
2679 AddIntVerticesInFaces( occgeo, *_ngMesh, nodeVec, internals );
2680 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2683 // Build viscous layers
2684 if ( _isViscousLayers2D )
2686 if ( !internals.hasInternalVertexInFace() ) {
2687 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2688 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2690 SMESH_ProxyMesh::Ptr viscousMesh;
2691 SMESH_MesherHelper helper( *_mesh );
2692 for ( int faceID = 1; faceID <= occgeo.fmap.Extent(); ++faceID )
2694 const TopoDS_Face& F = TopoDS::Face( occgeo.fmap( faceID ));
2695 viscousMesh = StdMeshers_ViscousLayers2D::Compute( *_mesh, F );
2698 // exclude from computation ng segments built on EDGEs of F
2699 for (int i = 1; i <= _ngMesh->GetNSeg(); i++)
2701 netgen::Segment & seg = _ngMesh->LineSegment(i);
2702 if (seg.si == faceID)
2705 // add new segments to _ngMesh instead of excluded ones
2706 helper.SetSubShape( F );
2708 StdMeshers_FaceSide::GetFaceWires( F, *_mesh, /*skipMediumNodes=*/true,
2709 error, viscousMesh );
2710 error = AddSegmentsToMesh( *_ngMesh, occgeo, wires, helper, nodeVec );
2712 if ( !error ) error = SMESH_ComputeError::New();
2714 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2717 // Let netgen compute 2D mesh
2718 startWith = netgen::MESHCONST_MESHSURFACE;
2719 endWith = _optimize ? netgen::MESHCONST_OPTSURFACE : netgen::MESHCONST_MESHSURFACE;
2724 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2726 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2728 if(netgen::multithread.terminate)
2731 comment << text (err);
2733 catch (Standard_Failure& ex)
2735 comment << text(ex);
2736 //err = 1; -- try to make volumes anyway
2738 catch (netgen::NgException exc)
2740 comment << text(exc);
2741 //err = 1; -- try to make volumes anyway
2746 doneTime += faceMeshingTime + ( _optimize ? faceOptimizTime : 0 );
2747 _ticTime = doneTime / _totalTime / _progressTic;
2749 // ---------------------
2750 // generate volume mesh
2751 // ---------------------
2752 // Fill _ngMesh with nodes and faces of computed 2D submeshes
2753 if ( !err && _isVolume && ( !meshedSM[ MeshDim_2D ].empty() || mparams.quad ))
2755 // load SMESH with computed segments and faces
2756 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2758 // compute pyramids on quadrangles
2759 SMESH_ProxyMesh::Ptr proxyMesh;
2760 if ( _mesh->NbQuadrangles() > 0 )
2761 for ( int iS = 1; iS <= occgeo.somap.Extent(); ++iS )
2763 StdMeshers_QuadToTriaAdaptor* Adaptor = new StdMeshers_QuadToTriaAdaptor;
2764 proxyMesh.reset( Adaptor );
2766 int nbPyrams = _mesh->NbPyramids();
2767 Adaptor->Compute( *_mesh, occgeo.somap(iS) );
2768 if ( nbPyrams != _mesh->NbPyramids() )
2770 list< SMESH_subMesh* > quadFaceSM;
2771 for (TopExp_Explorer face(occgeo.somap(iS), TopAbs_FACE); face.More(); face.Next())
2772 if ( Adaptor->GetProxySubMesh( face.Current() ))
2774 quadFaceSM.push_back( _mesh->GetSubMesh( face.Current() ));
2775 meshedSM[ MeshDim_2D ].remove( quadFaceSM.back() );
2777 FillNgMesh(occgeo, *_ngMesh, nodeVec, quadFaceSM, &quadHelper, proxyMesh);
2780 // fill _ngMesh with faces of sub-meshes
2781 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_2D ], &quadHelper));
2782 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2783 //toPython( _ngMesh, "/tmp/ngPython.py");
2785 if (!err && _isVolume)
2787 // Pass 3D simple parameters to NETGEN
2788 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
2789 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
2791 if ( double vol = simple3d->GetMaxElementVolume() ) {
2793 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
2794 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2797 // length from faces
2798 mparams.maxh = _ngMesh->AverageH();
2800 _ngMesh->SetGlobalH (mparams.maxh);
2801 mparams.grading = 0.4;
2803 _ngMesh->CalcLocalH(mparams.grading);
2805 _ngMesh->CalcLocalH();
2808 // Care of vertices internal in solids and internal faces (issue 0020676)
2809 if ( internals.hasInternalVertexInSolid() || internals.hasInternalFaces() )
2811 // store computed faces in SMESH in order not to create SMESH
2812 // faces for ng faces added here
2813 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2814 // add ng faces to solids with internal vertices
2815 AddIntVerticesInSolids( occgeo, *_ngMesh, nodeVec, internals );
2816 // duplicate mesh faces on internal faces
2817 FixIntFaces( occgeo, *_ngMesh, internals );
2818 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2820 // Let netgen compute 3D mesh
2821 startWith = endWith = netgen::MESHCONST_MESHVOLUME;
2826 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2828 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2830 if(netgen::multithread.terminate)
2833 if ( comment.empty() ) // do not overwrite a previos error
2834 comment << text(err);
2836 catch (Standard_Failure& ex)
2838 if ( comment.empty() ) // do not overwrite a previos error
2839 comment << text(ex);
2842 catch (netgen::NgException exc)
2844 if ( comment.empty() ) // do not overwrite a previos error
2845 comment << text(exc);
2848 _ticTime = ( doneTime += voluMeshingTime ) / _totalTime / _progressTic;
2850 // Let netgen optimize 3D mesh
2851 if ( !err && _optimize )
2853 startWith = endWith = netgen::MESHCONST_OPTVOLUME;
2858 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2860 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2862 if(netgen::multithread.terminate)
2865 if ( comment.empty() ) // do not overwrite a previos error
2866 comment << text(err);
2868 catch (Standard_Failure& ex)
2870 if ( comment.empty() ) // do not overwrite a previos error
2871 comment << text(ex);
2873 catch (netgen::NgException exc)
2875 if ( comment.empty() ) // do not overwrite a previos error
2876 comment << text(exc);
2880 if (!err && mparams.secondorder > 0)
2885 if ( !meshedSM[ MeshDim_1D ].empty() )
2887 // remove segments not attached to geometry (IPAL0052479)
2888 for (int i = 1; i <= _ngMesh->GetNSeg(); ++i)
2890 const netgen::Segment & seg = _ngMesh->LineSegment (i);
2891 if ( seg.epgeominfo[ 0 ].edgenr == 0 )
2892 _ngMesh->DeleteSegment( i );
2894 _ngMesh->Compress();
2896 // convert to quadratic
2897 netgen::OCCRefinementSurfaces ref (occgeo);
2898 ref.MakeSecondOrder (*_ngMesh);
2900 // care of elements already loaded to SMESH
2901 // if ( initState._nbSegments > 0 )
2902 // makeQuadratic( occgeo.emap, _mesh );
2903 // if ( initState._nbFaces > 0 )
2904 // makeQuadratic( occgeo.fmap, _mesh );
2906 catch (Standard_Failure& ex)
2908 if ( comment.empty() ) // do not overwrite a previos error
2909 comment << "Exception in netgen at passing to 2nd order ";
2911 catch (netgen::NgException exc)
2913 if ( comment.empty() ) // do not overwrite a previos error
2914 comment << exc.What();
2919 _ticTime = 0.98 / _progressTic;
2921 int nbNod = _ngMesh->GetNP();
2922 int nbSeg = _ngMesh->GetNSeg();
2923 int nbFac = _ngMesh->GetNSE();
2924 int nbVol = _ngMesh->GetNE();
2925 bool isOK = ( !err && (_isVolume ? (nbVol > 0) : (nbFac > 0)) );
2927 MESSAGE((err ? "Mesh Generation failure" : "End of Mesh Generation") <<
2928 ", nb nodes: " << nbNod <<
2929 ", nb segments: " << nbSeg <<
2930 ", nb faces: " << nbFac <<
2931 ", nb volumes: " << nbVol);
2933 // Feed back the SMESHDS with the generated Nodes and Elements
2934 if ( true /*isOK*/ ) // get whatever built
2936 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2938 if ( quadHelper.GetIsQuadratic() ) // remove free nodes
2939 for ( size_t i = 0; i < nodeVec.size(); ++i )
2940 if ( nodeVec[i] && nodeVec[i]->NbInverseElements() == 0 )
2941 _mesh->GetMeshDS()->RemoveFreeNode( nodeVec[i], 0, /*fromGroups=*/false );
2943 SMESH_ComputeErrorPtr readErr = ReadErrors(nodeVec);
2944 if ( readErr && !readErr->myBadElements.empty() )
2947 if ( !comment.empty() && !readErr->myComment.empty() ) comment += "\n";
2948 comment += readErr->myComment;
2950 if ( error->IsOK() && ( !isOK || comment.size() > 0 ))
2951 error->myName = COMPERR_ALGO_FAILED;
2952 if ( !comment.empty() )
2953 error->myComment = comment;
2955 // SetIsAlwaysComputed( true ) to empty sub-meshes, which
2956 // appear if the geometry contains coincident sub-shape due
2957 // to bool merge_solids = 1; in netgen/libsrc/occ/occgenmesh.cpp
2958 const int nbMaps = 2;
2959 const TopTools_IndexedMapOfShape* geoMaps[nbMaps] =
2960 { & occgeo.vmap, & occgeo.emap/*, & occgeo.fmap*/ };
2961 for ( int iMap = 0; iMap < nbMaps; ++iMap )
2962 for (int i = 1; i <= geoMaps[iMap]->Extent(); i++)
2963 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( geoMaps[iMap]->FindKey(i)))
2964 if ( !sm->IsMeshComputed() )
2965 sm->SetIsAlwaysComputed( true );
2967 // set bad compute error to subshapes of all failed sub-shapes
2968 if ( !error->IsOK() )
2970 bool pb2D = false, pb3D = false;
2971 for (int i = 1; i <= occgeo.fmap.Extent(); i++) {
2972 int status = occgeo.facemeshstatus[i-1];
2973 if (status == 1 ) continue;
2974 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.fmap( i ))) {
2975 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
2976 if ( !smError || smError->IsOK() ) {
2978 smError.reset( new SMESH_ComputeError( *error ));
2980 smError.reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED, "Ignored" ));
2981 if ( SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
2982 smError->myName = COMPERR_WARNING;
2984 pb2D = pb2D || smError->IsKO();
2987 if ( !pb2D ) // all faces are OK
2988 for (int i = 1; i <= occgeo.somap.Extent(); i++)
2989 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.somap( i )))
2991 bool smComputed = nbVol && !sm->IsEmpty();
2992 if ( smComputed && internals.hasInternalVertexInSolid( sm->GetId() ))
2994 int nbIntV = internals.getSolidsWithVertices().find( sm->GetId() )->second.size();
2995 SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
2996 smComputed = ( smDS->NbElements() > 0 || smDS->NbNodes() > nbIntV );
2998 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
2999 if ( !smComputed && ( !smError || smError->IsOK() ))
3001 smError.reset( new SMESH_ComputeError( *error ));
3002 if ( nbVol && SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
3004 smError->myName = COMPERR_WARNING;
3006 else if ( !smError->myBadElements.empty() ) // bad surface mesh
3008 if ( !hasBadElemOnSolid( smError->myBadElements, sm ))
3012 pb3D = pb3D || ( smError && smError->IsKO() );
3014 if ( !pb2D && !pb3D )
3015 err = 0; // no fatal errors, only warnings
3018 ngLib._isComputeOk = !err;
3023 //=============================================================================
3027 //=============================================================================
3028 bool NETGENPlugin_Mesher::Evaluate(MapShapeNbElems& aResMap)
3030 netgen::MeshingParameters& mparams = netgen::mparam;
3033 // -------------------------
3034 // Prepare OCC geometry
3035 // -------------------------
3036 netgen::OCCGeometry occgeo;
3037 list< SMESH_subMesh* > meshedSM[4]; // for 0-3 dimensions
3038 NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
3039 PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
3041 bool tooManyElems = false;
3042 const int hugeNb = std::numeric_limits<int>::max() / 100;
3047 // pass 1D simple parameters to NETGEN
3050 // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
3051 mparams.uselocalh = false;
3052 mparams.grading = 0.8; // not limitited size growth
3054 if ( _simpleHyp->GetNumberOfSegments() )
3056 mparams.maxh = occgeo.boundingbox.Diam();
3059 mparams.maxh = _simpleHyp->GetLocalLength();
3062 if ( mparams.maxh == 0.0 )
3063 mparams.maxh = occgeo.boundingbox.Diam();
3064 if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
3065 mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
3067 // let netgen create _ngMesh and calculate element size on not meshed shapes
3068 NETGENPlugin_NetgenLibWrapper ngLib;
3069 netgen::Mesh *ngMesh = NULL;
3073 int startWith = netgen::MESHCONST_ANALYSE;
3074 int endWith = netgen::MESHCONST_MESHEDGES;
3076 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, mparams, startWith, endWith);
3078 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, startWith, endWith, optstr);
3081 if(netgen::multithread.terminate)
3084 ngLib.setMesh(( Ng_Mesh*) ngMesh );
3086 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( _shape ))
3087 sm->GetComputeError().reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED ));
3092 // Pass 1D simple parameters to NETGEN
3093 // --------------------------------
3094 int nbSeg = _simpleHyp->GetNumberOfSegments();
3095 double segSize = _simpleHyp->GetLocalLength();
3096 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
3098 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
3100 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
3101 setLocalSize( e, segSize, *ngMesh );
3104 else // if ( ! _simpleHyp )
3106 // Local size on vertices and edges
3107 // --------------------------------
3108 for(std::map<int,double>::const_iterator it=EdgeId2LocalSize.begin(); it!=EdgeId2LocalSize.end(); it++)
3110 int key = (*it).first;
3111 double hi = (*it).second;
3112 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3113 const TopoDS_Edge& e = TopoDS::Edge(shape);
3114 setLocalSize( e, hi, *ngMesh );
3116 for(std::map<int,double>::const_iterator it=VertexId2LocalSize.begin(); it!=VertexId2LocalSize.end(); it++)
3118 int key = (*it).first;
3119 double hi = (*it).second;
3120 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3121 const TopoDS_Vertex& v = TopoDS::Vertex(shape);
3122 gp_Pnt p = BRep_Tool::Pnt(v);
3123 NETGENPlugin_Mesher::RestrictLocalSize( *ngMesh, p.XYZ(), hi );
3125 for(map<int,double>::const_iterator it=FaceId2LocalSize.begin();
3126 it!=FaceId2LocalSize.end(); it++)
3128 int key = (*it).first;
3129 double val = (*it).second;
3130 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3131 int faceNgID = occgeo.fmap.FindIndex(shape);
3132 occgeo.SetFaceMaxH(faceNgID, val);
3133 for ( TopExp_Explorer edgeExp( shape, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
3134 setLocalSize( TopoDS::Edge( edgeExp.Current() ), val, *ngMesh );
3137 // calculate total nb of segments and length of edges
3138 double fullLen = 0.0;
3140 int entity = mparams.secondorder > 0 ? SMDSEntity_Quad_Edge : SMDSEntity_Edge;
3141 TopTools_DataMapOfShapeInteger Edge2NbSeg;
3142 for (TopExp_Explorer exp(_shape, TopAbs_EDGE); exp.More(); exp.Next())
3144 TopoDS_Edge E = TopoDS::Edge( exp.Current() );
3145 if( !Edge2NbSeg.Bind(E,0) )
3148 double aLen = SMESH_Algo::EdgeLength(E);
3151 vector<int>& aVec = aResMap[_mesh->GetSubMesh(E)];
3153 aVec.resize( SMDSEntity_Last, 0);
3155 fullNbSeg += aVec[ entity ];
3158 // store nb of segments computed by Netgen
3159 NCollection_Map<Link> linkMap;
3160 for (int i = 1; i <= ngMesh->GetNSeg(); ++i )
3162 const netgen::Segment& seg = ngMesh->LineSegment(i);
3163 Link link(seg[0], seg[1]);
3164 if ( !linkMap.Add( link )) continue;
3165 int aGeomEdgeInd = seg.epgeominfo[0].edgenr;
3166 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
3168 vector<int>& aVec = aResMap[_mesh->GetSubMesh(occgeo.emap(aGeomEdgeInd))];
3172 // store nb of nodes on edges computed by Netgen
3173 TopTools_DataMapIteratorOfDataMapOfShapeInteger Edge2NbSegIt(Edge2NbSeg);
3174 for (; Edge2NbSegIt.More(); Edge2NbSegIt.Next())
3176 vector<int>& aVec = aResMap[_mesh->GetSubMesh(Edge2NbSegIt.Key())];
3177 if ( aVec[ entity ] > 1 && aVec[ SMDSEntity_Node ] == 0 )
3178 aVec[SMDSEntity_Node] = mparams.secondorder > 0 ? 2*aVec[ entity ]-1 : aVec[ entity ]-1;
3180 fullNbSeg += aVec[ entity ];
3181 Edge2NbSeg( Edge2NbSegIt.Key() ) = aVec[ entity ];
3183 if ( fullNbSeg == 0 )
3190 if ( double area = _simpleHyp->GetMaxElementArea() ) {
3192 mparams.maxh = sqrt(2. * area/sqrt(3.0));
3193 mparams.grading = 0.4; // moderate size growth
3196 // length from edges
3197 mparams.maxh = fullLen/fullNbSeg;
3198 mparams.grading = 0.2; // slow size growth
3201 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3202 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3204 for (TopExp_Explorer exp(_shape, TopAbs_FACE); exp.More(); exp.Next())
3206 TopoDS_Face F = TopoDS::Face( exp.Current() );
3207 SMESH_subMesh *sm = _mesh->GetSubMesh(F);
3209 BRepGProp::SurfaceProperties(F,G);
3210 double anArea = G.Mass();
3211 tooManyElems = tooManyElems || ( anArea/hugeNb > mparams.maxh*mparams.maxh );
3213 if ( !tooManyElems )
3215 TopTools_MapOfShape egdes;
3216 for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next())
3217 if ( egdes.Add( exp1.Current() ))
3218 nb1d += Edge2NbSeg.Find(exp1.Current());
3220 int nbFaces = tooManyElems ? hugeNb : int( 4*anArea / (mparams.maxh*mparams.maxh*sqrt(3.)));
3221 int nbNodes = tooManyElems ? hugeNb : (( nbFaces*3 - (nb1d-1)*2 ) / 6 + 1 );
3223 vector<int> aVec(SMDSEntity_Last, 0);
3224 if( mparams.secondorder > 0 ) {
3225 int nb1d_in = (nbFaces*3 - nb1d) / 2;
3226 aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
3227 aVec[SMDSEntity_Quad_Triangle] = nbFaces;
3230 aVec[SMDSEntity_Node] = Max ( nbNodes, 0 );
3231 aVec[SMDSEntity_Triangle] = nbFaces;
3233 aResMap[sm].swap(aVec);
3240 // pass 3D simple parameters to NETGEN
3241 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
3242 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
3244 if ( double vol = simple3d->GetMaxElementVolume() ) {
3246 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
3247 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3250 // using previous length from faces
3252 mparams.grading = 0.4;
3253 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3256 BRepGProp::VolumeProperties(_shape,G);
3257 double aVolume = G.Mass();
3258 double tetrVol = 0.1179*mparams.maxh*mparams.maxh*mparams.maxh;
3259 tooManyElems = tooManyElems || ( aVolume/hugeNb > tetrVol );
3260 int nbVols = tooManyElems ? hugeNb : int(aVolume/tetrVol);
3261 int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
3262 vector<int> aVec(SMDSEntity_Last, 0 );
3263 if ( tooManyElems ) // avoid FPE
3265 aVec[SMDSEntity_Node] = hugeNb;
3266 aVec[ mparams.secondorder > 0 ? SMDSEntity_Quad_Tetra : SMDSEntity_Tetra] = hugeNb;
3270 if( mparams.secondorder > 0 ) {
3271 aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
3272 aVec[SMDSEntity_Quad_Tetra] = nbVols;
3275 aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
3276 aVec[SMDSEntity_Tetra] = nbVols;
3279 SMESH_subMesh *sm = _mesh->GetSubMesh(_shape);
3280 aResMap[sm].swap(aVec);
3286 double NETGENPlugin_Mesher::GetProgress(const SMESH_Algo* holder,
3287 const int * algoProgressTic,
3288 const double * algoProgress) const
3290 ((int&) _progressTic ) = *algoProgressTic + 1;
3292 if ( !_occgeom ) return 0;
3294 double progress = -1;
3297 if ( _ticTime < 0 && netgen::multithread.task[0] == 'O'/*Optimizing surface*/ )
3299 ((double&) _ticTime ) = edgeFaceMeshingTime / _totalTime / _progressTic;
3301 else if ( !_optimize /*&& _occgeom->fmap.Extent() > 1*/ )
3303 int doneShapeIndex = -1;
3304 while ( doneShapeIndex+1 < _occgeom->facemeshstatus.Size() &&
3305 _occgeom->facemeshstatus[ doneShapeIndex+1 ])
3307 if ( doneShapeIndex+1 != _curShapeIndex )
3309 ((int&) _curShapeIndex) = doneShapeIndex+1;
3310 double doneShapeRate = _curShapeIndex / double( _occgeom->fmap.Extent() );
3311 double doneTime = edgeMeshingTime + doneShapeRate * faceMeshingTime;
3312 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3313 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3314 // << " " << doneTime / _totalTime / _progressTic << endl;
3318 else if ( !_optimize && _occgeom->somap.Extent() > 1 )
3320 int curShapeIndex = _curShapeIndex;
3321 if ( _ngMesh->GetNE() > 0 )
3323 netgen::Element el = (*_ngMesh)[netgen::ElementIndex( _ngMesh->GetNE()-1 )];
3324 curShapeIndex = el.GetIndex();
3326 if ( curShapeIndex != _curShapeIndex )
3328 ((int&) _curShapeIndex) = curShapeIndex;
3329 double doneShapeRate = _curShapeIndex / double( _occgeom->somap.Extent() );
3330 double doneTime = edgeFaceMeshingTime + doneShapeRate * voluMeshingTime;
3331 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3332 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3333 // << " " << doneTime / _totalTime / _progressTic << endl;
3337 progress = Max( *algoProgressTic * _ticTime, *algoProgress );
3340 ((int&) *algoProgressTic )++;
3341 ((double&) *algoProgress) = progress;
3343 //cout << progress << " " << *algoProgressTic << " " << netgen::multithread.task << " "<< _ticTime << endl;
3345 return Min( progress, 0.99 );
3348 //================================================================================
3350 * \brief Remove "test.out" and "problemfaces" files in current directory
3352 //================================================================================
3354 void NETGENPlugin_Mesher::RemoveTmpFiles()
3356 bool rm = SMESH_File("test.out").remove() ;
3358 if (rm && netgen::testout)
3360 delete netgen::testout;
3361 netgen::testout = 0;
3364 SMESH_File("problemfaces").remove();
3365 SMESH_File("occmesh.rep").remove();
3368 //================================================================================
3370 * \brief Read mesh entities preventing successful computation from "test.out" file
3372 //================================================================================
3374 SMESH_ComputeErrorPtr
3375 NETGENPlugin_Mesher::ReadErrors(const vector<const SMDS_MeshNode* >& nodeVec)
3377 SMESH_ComputeErrorPtr err = SMESH_ComputeError::New
3378 (COMPERR_BAD_INPUT_MESH, "Some edges multiple times in surface mesh");
3379 SMESH_File file("test.out");
3381 vector<int> three1(3), three2(3);
3382 const char* badEdgeStr = " multiple times in surface mesh";
3383 const int badEdgeStrLen = strlen( badEdgeStr );
3385 while( !file.eof() )
3387 if ( strncmp( file, "Edge ", 5 ) == 0 &&
3388 file.getInts( two ) &&
3389 strncmp( file, badEdgeStr, badEdgeStrLen ) == 0 &&
3390 two[0] < nodeVec.size() && two[1] < nodeVec.size())
3392 err->myBadElements.push_back( new SMDS_LinearEdge( nodeVec[ two[0]], nodeVec[ two[1]] ));
3393 file += badEdgeStrLen;
3395 else if ( strncmp( file, "Intersecting: ", 14 ) == 0 )
3398 // openelement 18 with open element 126
3402 const char* pos = file;
3403 bool ok = ( strncmp( file, "openelement ", 12 ) == 0 );
3404 ok = ok && file.getInts( two );
3405 ok = ok && file.getInts( three1 );
3406 ok = ok && file.getInts( three2 );
3407 for ( int i = 0; ok && i < 3; ++i )
3408 ok = ( three1[i] < nodeVec.size() && nodeVec[ three1[i]]);
3409 for ( int i = 0; ok && i < 3; ++i )
3410 ok = ( three2[i] < nodeVec.size() && nodeVec[ three2[i]]);
3413 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three1[0]],
3414 nodeVec[ three1[1]],
3415 nodeVec[ three1[2]]));
3416 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three2[0]],
3417 nodeVec[ three2[1]],
3418 nodeVec[ three2[2]]));
3419 err->myComment = "Intersecting triangles";
3433 size_t nbBadElems = err->myBadElements.size();
3440 //================================================================================
3442 * \brief Write a python script creating an equivalent SALOME mesh.
3443 * This is useful to see what mesh is passed as input for the next step of mesh
3444 * generation (of mesh of higher dimension)
3446 //================================================================================
3448 void NETGENPlugin_Mesher::toPython( const netgen::Mesh* ngMesh,
3449 const std::string& pyFile)
3451 ofstream outfile(pyFile.c_str(), ios::out);
3452 if ( !outfile ) return;
3454 outfile << "import SMESH" << endl
3455 << "from salome.smesh import smeshBuilder" << endl
3456 << "smesh = smeshBuilder.New(salome.myStudy)" << endl
3457 << "mesh = smesh.Mesh()" << endl << endl;
3459 using namespace netgen;
3461 for (pi = PointIndex::BASE;
3462 pi < ngMesh->GetNP()+PointIndex::BASE; pi++)
3464 outfile << "mesh.AddNode( ";
3465 outfile << (*ngMesh)[pi](0) << ", ";
3466 outfile << (*ngMesh)[pi](1) << ", ";
3467 outfile << (*ngMesh)[pi](2) << ") ## "<< pi << endl;
3470 int nbDom = ngMesh->GetNDomains();
3471 for ( int i = 0; i < nbDom; ++i )
3472 outfile<< "grp" << i+1 << " = mesh.CreateEmptyGroup( SMESH.FACE, 'domain"<< i+1 << "')"<< endl;
3474 SurfaceElementIndex sei;
3475 for (sei = 0; sei < ngMesh->GetNSE(); sei++)
3477 outfile << "mesh.AddFace([ ";
3478 Element2d sel = (*ngMesh)[sei];
3479 for (int j = 0; j < sel.GetNP(); j++)
3480 outfile << sel[j] << ( j+1 < sel.GetNP() ? ", " : " ])");
3481 if ( sel.IsDeleted() ) outfile << " ## IsDeleted ";
3484 if ((*ngMesh)[sei].GetIndex())
3486 if ( int dom1 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainIn())
3487 outfile << "grp"<< dom1 <<".Add([ " << (int)sei+1 << " ])" << endl;
3488 if ( int dom2 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainOut())
3489 outfile << "grp"<< dom2 <<".Add([ " << (int)sei+1 << " ])" << endl;
3493 for (ElementIndex ei = 0; ei < ngMesh->GetNE(); ei++)
3495 Element el = (*ngMesh)[ei];
3496 outfile << "mesh.AddVolume([ ";
3497 for (int j = 0; j < el.GetNP(); j++)
3498 outfile << el[j] << ( j+1 < el.GetNP() ? ", " : " ])");
3502 for (int i = 1; i <= ngMesh->GetNSeg(); i++)
3504 const Segment & seg = ngMesh->LineSegment (i);
3505 outfile << "mesh.AddEdge([ "
3507 << seg[1] << " ])" << endl;
3509 cout << "Write " << pyFile << endl;
3512 //================================================================================
3514 * \brief Constructor of NETGENPlugin_ngMeshInfo
3516 //================================================================================
3518 NETGENPlugin_ngMeshInfo::NETGENPlugin_ngMeshInfo( netgen::Mesh* ngMesh):
3523 _nbNodes = ngMesh->GetNP();
3524 _nbSegments = ngMesh->GetNSeg();
3525 _nbFaces = ngMesh->GetNSE();
3526 _nbVolumes = ngMesh->GetNE();
3530 _nbNodes = _nbSegments = _nbFaces = _nbVolumes = 0;
3534 //================================================================================
3536 * \brief Copy LocalH member from one netgen mesh to another
3538 //================================================================================
3540 void NETGENPlugin_ngMeshInfo::transferLocalH( netgen::Mesh* fromMesh,
3541 netgen::Mesh* toMesh )
3543 if ( !fromMesh->LocalHFunctionGenerated() ) return;
3544 if ( !toMesh->LocalHFunctionGenerated() )
3546 toMesh->CalcLocalH(netgen::mparam.grading);
3548 toMesh->CalcLocalH();
3551 const size_t size = sizeof( netgen::LocalH );
3552 _copyOfLocalH = new char[ size ];
3553 memcpy( (void*)_copyOfLocalH, (void*)&toMesh->LocalHFunction(), size );
3554 memcpy( (void*)&toMesh->LocalHFunction(), (void*)&fromMesh->LocalHFunction(), size );
3557 //================================================================================
3559 * \brief Restore LocalH member of a netgen mesh
3561 //================================================================================
3563 void NETGENPlugin_ngMeshInfo::restoreLocalH( netgen::Mesh* toMesh )
3565 if ( _copyOfLocalH )
3567 const size_t size = sizeof( netgen::LocalH );
3568 memcpy( (void*)&toMesh->LocalHFunction(), (void*)_copyOfLocalH, size );
3569 delete [] _copyOfLocalH;
3574 //================================================================================
3576 * \brief Find "internal" sub-shapes
3578 //================================================================================
3580 NETGENPlugin_Internals::NETGENPlugin_Internals( SMESH_Mesh& mesh,
3581 const TopoDS_Shape& shape,
3583 : _mesh( mesh ), _is3D( is3D )
3585 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3587 TopExp_Explorer f,e;
3588 for ( f.Init( shape, TopAbs_FACE ); f.More(); f.Next() )
3590 int faceID = meshDS->ShapeToIndex( f.Current() );
3592 // find not computed internal edges
3594 for ( e.Init( f.Current().Oriented(TopAbs_FORWARD), TopAbs_EDGE ); e.More(); e.Next() )
3595 if ( e.Current().Orientation() == TopAbs_INTERNAL )
3597 SMESH_subMesh* eSM = mesh.GetSubMesh( e.Current() );
3598 if ( eSM->IsEmpty() )
3600 _e2face.insert( make_pair( eSM->GetId(), faceID ));
3601 for ( TopoDS_Iterator v(e.Current()); v.More(); v.Next() )
3602 _e2face.insert( make_pair( meshDS->ShapeToIndex( v.Value() ), faceID ));
3606 // find internal vertices in a face
3607 set<int> intVV; // issue 0020850 where same vertex is twice in a face
3608 for ( TopoDS_Iterator fSub( f.Current() ); fSub.More(); fSub.Next())
3609 if ( fSub.Value().ShapeType() == TopAbs_VERTEX )
3611 int vID = meshDS->ShapeToIndex( fSub.Value() );
3612 if ( intVV.insert( vID ).second )
3613 _f2v[ faceID ].push_back( vID );
3618 // find internal faces and their subshapes where nodes are to be doubled
3619 // to make a crack with non-sewed borders
3621 if ( f.Current().Orientation() == TopAbs_INTERNAL )
3623 _intShapes.insert( meshDS->ShapeToIndex( f.Current() ));
3626 list< TopoDS_Shape > edges;
3627 for ( e.Init( f.Current(), TopAbs_EDGE ); e.More(); e.Next())
3628 if ( SMESH_MesherHelper::NbAncestors( e.Current(), mesh, TopAbs_FACE ) > 1 )
3630 _intShapes.insert( meshDS->ShapeToIndex( e.Current() ));
3631 edges.push_back( e.Current() );
3632 // find border faces
3633 PShapeIteratorPtr fIt =
3634 SMESH_MesherHelper::GetAncestors( edges.back(),mesh,TopAbs_FACE );
3635 while ( const TopoDS_Shape* pFace = fIt->next() )
3636 if ( !pFace->IsSame( f.Current() ))
3637 _borderFaces.insert( meshDS->ShapeToIndex( *pFace ));
3640 // we consider vertex internal if it is shared by more than one internal edge
3641 list< TopoDS_Shape >::iterator edge = edges.begin();
3642 for ( ; edge != edges.end(); ++edge )
3643 for ( TopoDS_Iterator v( *edge ); v.More(); v.Next() )
3645 set<int> internalEdges;
3646 PShapeIteratorPtr eIt =
3647 SMESH_MesherHelper::GetAncestors( v.Value(),mesh,TopAbs_EDGE );
3648 while ( const TopoDS_Shape* pEdge = eIt->next() )
3650 int edgeID = meshDS->ShapeToIndex( *pEdge );
3651 if ( isInternalShape( edgeID ))
3652 internalEdges.insert( edgeID );
3654 if ( internalEdges.size() > 1 )
3655 _intShapes.insert( meshDS->ShapeToIndex( v.Value() ));
3659 } // loop on geom faces
3661 // find vertices internal in solids
3664 for ( TopExp_Explorer so(shape, TopAbs_SOLID); so.More(); so.Next())
3666 int soID = meshDS->ShapeToIndex( so.Current() );
3667 for ( TopoDS_Iterator soSub( so.Current() ); soSub.More(); soSub.Next())
3668 if ( soSub.Value().ShapeType() == TopAbs_VERTEX )
3669 _s2v[ soID ].push_back( meshDS->ShapeToIndex( soSub.Value() ));
3674 //================================================================================
3676 * \brief Find mesh faces on non-internal geom faces sharing internal edge
3677 * some nodes of which are to be doubled to make the second border of the "crack"
3679 //================================================================================
3681 void NETGENPlugin_Internals::findBorderElements( TIDSortedElemSet & borderElems )
3683 if ( _intShapes.empty() ) return;
3685 SMESH_Mesh& mesh = const_cast<SMESH_Mesh&>(_mesh);
3686 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3688 // loop on internal geom edges
3689 set<int>::const_iterator intShapeId = _intShapes.begin();
3690 for ( ; intShapeId != _intShapes.end(); ++intShapeId )
3692 const TopoDS_Shape& s = meshDS->IndexToShape( *intShapeId );
3693 if ( s.ShapeType() != TopAbs_EDGE ) continue;
3695 // get internal and non-internal geom faces sharing the internal edge <s>
3697 set<int>::iterator bordFace = _borderFaces.end();
3698 PShapeIteratorPtr faces = SMESH_MesherHelper::GetAncestors( s, _mesh, TopAbs_FACE );
3699 while ( const TopoDS_Shape* pFace = faces->next() )
3701 int faceID = meshDS->ShapeToIndex( *pFace );
3702 if ( isInternalShape( faceID ))
3705 bordFace = _borderFaces.insert( faceID ).first;
3707 if ( bordFace == _borderFaces.end() || !intFace ) continue;
3709 // get all links of mesh faces on internal geom face sharing nodes on edge <s>
3710 set< SMESH_OrientedLink > links; //!< links of faces on internal geom face
3711 list<const SMDS_MeshElement*> suspectFaces[2]; //!< mesh faces on border geom faces
3712 int nbSuspectFaces = 0;
3713 SMESHDS_SubMesh* intFaceSM = meshDS->MeshElements( intFace );
3714 if ( !intFaceSM || intFaceSM->NbElements() == 0 ) continue;
3715 SMESH_subMeshIteratorPtr smIt = mesh.GetSubMesh( s )->getDependsOnIterator(true,true);
3716 while ( smIt->more() )
3718 SMESHDS_SubMesh* sm = smIt->next()->GetSubMeshDS();
3719 if ( !sm ) continue;
3720 SMDS_NodeIteratorPtr nIt = sm->GetNodes();
3721 while ( nIt->more() )
3723 const SMDS_MeshNode* nOnEdge = nIt->next();
3724 SMDS_ElemIteratorPtr fIt = nOnEdge->GetInverseElementIterator(SMDSAbs_Face);
3725 while ( fIt->more() )
3727 const SMDS_MeshElement* f = fIt->next();
3728 const int nbNodes = f->NbCornerNodes();
3729 if ( intFaceSM->Contains( f ))
3731 for ( int i = 0; i < nbNodes; ++i )
3732 links.insert( SMESH_OrientedLink( f->GetNode(i), f->GetNode((i+1)%nbNodes)));
3737 for ( int i = 0; i < nbNodes; ++i )
3738 nbDblNodes += isInternalShape( f->GetNode(i)->getshapeId() );
3740 suspectFaces[ nbDblNodes < 2 ].push_back( f );
3746 // suspectFaces[0] having link with same orientation as mesh faces on
3747 // the internal geom face are <borderElems>. suspectFaces[1] have
3748 // only one node on edge <s>, we decide on them later (at the 2nd loop)
3749 // by links of <borderElems> found at the 1st and 2nd loops
3750 set< SMESH_OrientedLink > borderLinks;
3751 for ( int isPostponed = 0; isPostponed < 2; ++isPostponed )
3753 list<const SMDS_MeshElement*>::iterator fIt = suspectFaces[isPostponed].begin();
3754 for ( int nbF = 0; fIt != suspectFaces[isPostponed].end(); ++fIt, ++nbF )
3756 const SMDS_MeshElement* f = *fIt;
3757 bool isBorder = false, linkFound = false, borderLinkFound = false;
3758 list< SMESH_OrientedLink > faceLinks;
3759 int nbNodes = f->NbCornerNodes();
3760 for ( int i = 0; i < nbNodes; ++i )
3762 SMESH_OrientedLink link( f->GetNode(i), f->GetNode((i+1)%nbNodes));
3763 faceLinks.push_back( link );
3766 set< SMESH_OrientedLink >::iterator foundLink = links.find( link );
3767 if ( foundLink != links.end() )
3770 isBorder = ( foundLink->_reversed == link._reversed );
3771 if ( !isBorder && !isPostponed ) break;
3772 faceLinks.pop_back();
3774 else if ( isPostponed && !borderLinkFound )
3776 foundLink = borderLinks.find( link );
3777 if ( foundLink != borderLinks.end() )
3779 borderLinkFound = true;
3780 isBorder = ( foundLink->_reversed != link._reversed );
3787 borderElems.insert( f );
3788 borderLinks.insert( faceLinks.begin(), faceLinks.end() );
3790 else if ( !linkFound && !borderLinkFound )
3792 suspectFaces[1].push_back( f );
3793 if ( nbF > 2 * nbSuspectFaces )
3794 break; // dead loop protection
3801 //================================================================================
3803 * \brief put internal shapes in maps and fill in submeshes to precompute
3805 //================================================================================
3807 void NETGENPlugin_Internals::getInternalEdges( TopTools_IndexedMapOfShape& fmap,
3808 TopTools_IndexedMapOfShape& emap,
3809 TopTools_IndexedMapOfShape& vmap,
3810 list< SMESH_subMesh* > smToPrecompute[])
3812 if ( !hasInternalEdges() ) return;
3813 map<int,int>::const_iterator ev_face = _e2face.begin();
3814 for ( ; ev_face != _e2face.end(); ++ev_face )
3816 const TopoDS_Shape& ev = _mesh.GetMeshDS()->IndexToShape( ev_face->first );
3817 const TopoDS_Shape& face = _mesh.GetMeshDS()->IndexToShape( ev_face->second );
3819 ( ev.ShapeType() == TopAbs_EDGE ? emap : vmap ).Add( ev );
3821 //cout<<"INTERNAL EDGE or VERTEX "<<ev_face->first<<" on face "<<ev_face->second<<endl;
3823 smToPrecompute[ MeshDim_1D ].push_back( _mesh.GetSubMeshContaining( ev_face->first ));
3827 //================================================================================
3829 * \brief return shapes and submeshes to be meshed and already meshed boundary submeshes
3831 //================================================================================
3833 void NETGENPlugin_Internals::getInternalFaces( TopTools_IndexedMapOfShape& fmap,
3834 TopTools_IndexedMapOfShape& emap,
3835 list< SMESH_subMesh* >& intFaceSM,
3836 list< SMESH_subMesh* >& boundarySM)
3838 if ( !hasInternalFaces() ) return;
3840 // <fmap> and <emap> are for not yet meshed shapes
3841 // <intFaceSM> is for submeshes of faces
3842 // <boundarySM> is for meshed edges and vertices
3847 set<int> shapeIDs ( _intShapes );
3848 if ( !_borderFaces.empty() )
3849 shapeIDs.insert( _borderFaces.begin(), _borderFaces.end() );
3851 set<int>::const_iterator intS = shapeIDs.begin();
3852 for ( ; intS != shapeIDs.end(); ++intS )
3854 SMESH_subMesh* sm = _mesh.GetSubMeshContaining( *intS );
3856 if ( sm->GetSubShape().ShapeType() != TopAbs_FACE ) continue;
3858 intFaceSM.push_back( sm );
3860 // add submeshes of not computed internal faces
3861 if ( !sm->IsEmpty() ) continue;
3863 SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(true,true);
3864 while ( smIt->more() )
3867 const TopoDS_Shape& s = sm->GetSubShape();
3869 if ( sm->IsEmpty() )
3872 switch ( s.ShapeType() ) {
3873 case TopAbs_FACE: fmap.Add ( s ); break;
3874 case TopAbs_EDGE: emap.Add ( s ); break;
3880 if ( s.ShapeType() != TopAbs_FACE )
3881 boundarySM.push_back( sm );
3887 //================================================================================
3889 * \brief Return true if given shape is to be precomputed in order to be correctly
3890 * added to netgen mesh
3892 //================================================================================
3894 bool NETGENPlugin_Internals::isShapeToPrecompute(const TopoDS_Shape& s)
3896 int shapeID = _mesh.GetMeshDS()->ShapeToIndex( s );
3897 switch ( s.ShapeType() ) {
3898 case TopAbs_FACE : break; //return isInternalShape( shapeID ) || isBorderFace( shapeID );
3899 case TopAbs_EDGE : return isInternalEdge( shapeID );
3900 case TopAbs_VERTEX: break;
3906 //================================================================================
3908 * \brief Return SMESH
3910 //================================================================================
3912 SMESH_Mesh& NETGENPlugin_Internals::getMesh() const
3914 return const_cast<SMESH_Mesh&>( _mesh );
3917 //================================================================================
3919 * \brief Initialize netgen library
3921 //================================================================================
3923 NETGENPlugin_NetgenLibWrapper::NETGENPlugin_NetgenLibWrapper()
3927 _isComputeOk = false;
3929 if ( !getenv( "KEEP_NETGEN_OUTPUT" ))
3931 // redirect all netgen output (mycout,myerr,cout) to _outputFileName
3932 _outputFileName = getOutputFileName();
3933 netgen::mycout = new ofstream ( _outputFileName.c_str() );
3934 netgen::myerr = netgen::mycout;
3935 _coutBuffer = std::cout.rdbuf();
3937 cout << "NOTE: netgen output is redirected to file " << _outputFileName << endl;
3939 std::cout.rdbuf( netgen::mycout->rdbuf() );
3943 _ngMesh = Ng_NewMesh();
3946 //================================================================================
3948 * \brief Finish using netgen library
3950 //================================================================================
3952 NETGENPlugin_NetgenLibWrapper::~NETGENPlugin_NetgenLibWrapper()
3954 Ng_DeleteMesh( _ngMesh );
3956 NETGENPlugin_Mesher::RemoveTmpFiles();
3958 std::cout.rdbuf( _coutBuffer );
3965 //================================================================================
3967 * \brief Set netgen mesh to delete at destruction
3969 //================================================================================
3971 void NETGENPlugin_NetgenLibWrapper::setMesh( Ng_Mesh* mesh )
3974 Ng_DeleteMesh( _ngMesh );
3978 //================================================================================
3980 * \brief Return a unique file name
3982 //================================================================================
3984 std::string NETGENPlugin_NetgenLibWrapper::getOutputFileName()
3986 std::string aTmpDir = SALOMEDS_Tool::GetTmpDir();
3988 TCollection_AsciiString aGenericName = (char*)aTmpDir.c_str();
3989 aGenericName += "NETGEN_";
3991 aGenericName += getpid();
3993 aGenericName += _getpid();
3995 aGenericName += "_";
3996 aGenericName += Abs((Standard_Integer)(long) aGenericName.ToCString());
3997 aGenericName += ".out";
3999 return aGenericName.ToCString();
4002 //================================================================================
4004 * \brief Remove file with netgen output
4006 //================================================================================
4008 void NETGENPlugin_NetgenLibWrapper::removeOutputFile()
4010 if ( !_outputFileName.empty() )
4012 if ( netgen::mycout )
4014 delete netgen::mycout;
4018 string tmpDir = SALOMEDS_Tool::GetDirFromPath ( _outputFileName );
4019 string aFileName = SALOMEDS_Tool::GetNameFromPath( _outputFileName ) + ".out";
4020 SALOMEDS::ListOfFileNames_var aFiles = new SALOMEDS::ListOfFileNames;
4022 aFiles[0] = aFileName.c_str();
4024 SALOMEDS_Tool::RemoveTemporaryFiles( tmpDir.c_str(), aFiles.in(), true );