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>
60 #include <TopExp_Explorer.hxx>
61 #include <TopTools_DataMapIteratorOfDataMapOfShapeInteger.hxx>
62 #include <TopTools_DataMapIteratorOfDataMapOfShapeShape.hxx>
63 #include <TopTools_DataMapOfShapeInteger.hxx>
64 #include <TopTools_DataMapOfShapeShape.hxx>
65 #include <TopTools_MapOfShape.hxx>
67 #include <OSD_File.hxx>
68 #include <OSD_Path.hxx>
70 // Netgen include files
74 #include <occgeom.hpp>
75 #include <meshing.hpp>
76 //#include <ngexception.hpp>
79 extern int OCCGenerateMesh (OCCGeometry&, Mesh*&, MeshingParameters&, int, int);
81 extern int OCCGenerateMesh (OCCGeometry&, Mesh*&, int, int, char*);
83 //extern void OCCSetLocalMeshSize(OCCGeometry & geom, Mesh & mesh);
84 extern MeshingParameters mparam;
85 extern volatile multithreadt multithread;
86 extern bool merge_solids;
95 using namespace nglib;
99 #define nodeVec_ACCESS(index) ((SMDS_MeshNode*) nodeVec.at((index)))
101 #define nodeVec_ACCESS(index) ((SMDS_MeshNode*) nodeVec[index])
104 #define NGPOINT_COORDS(p) p(0),p(1),p(2)
107 // dump elements added to ng mesh
108 //#define DUMP_SEGMENTS
109 //#define DUMP_TRIANGLES
110 //#define DUMP_TRIANGLES_SCRIPT "/tmp/trias.py" //!< debug AddIntVerticesInSolids()
113 TopTools_IndexedMapOfShape ShapesWithLocalSize;
114 std::map<int,double> VertexId2LocalSize;
115 std::map<int,double> EdgeId2LocalSize;
116 std::map<int,double> FaceId2LocalSize;
118 //=============================================================================
122 //=============================================================================
124 NETGENPlugin_Mesher::NETGENPlugin_Mesher (SMESH_Mesh* mesh,
125 const TopoDS_Shape& aShape,
131 _fineness(NETGENPlugin_Hypothesis::GetDefaultFineness()),
132 _isViscousLayers2D(false),
141 SetDefaultParameters();
142 ShapesWithLocalSize.Clear();
143 VertexId2LocalSize.clear();
144 EdgeId2LocalSize.clear();
145 FaceId2LocalSize.clear();
148 //================================================================================
152 //================================================================================
154 NETGENPlugin_Mesher::~NETGENPlugin_Mesher()
162 //================================================================================
164 * Set pointer to NETGENPlugin_Mesher* field of the holder, that will be
165 * nullified at destruction of this
167 //================================================================================
169 void NETGENPlugin_Mesher::SetSelfPointer( NETGENPlugin_Mesher ** ptr )
180 //================================================================================
182 * \brief Initialize global NETGEN parameters with default values
184 //================================================================================
186 void NETGENPlugin_Mesher::SetDefaultParameters()
188 netgen::MeshingParameters& mparams = netgen::mparam;
189 // maximal mesh edge size
190 mparams.maxh = 0;//NETGENPlugin_Hypothesis::GetDefaultMaxSize();
192 // minimal number of segments per edge
193 mparams.segmentsperedge = NETGENPlugin_Hypothesis::GetDefaultNbSegPerEdge();
194 // rate of growth of size between elements
195 mparams.grading = NETGENPlugin_Hypothesis::GetDefaultGrowthRate();
196 // safety factor for curvatures (elements per radius)
197 mparams.curvaturesafety = NETGENPlugin_Hypothesis::GetDefaultNbSegPerRadius();
198 // create elements of second order
199 mparams.secondorder = NETGENPlugin_Hypothesis::GetDefaultSecondOrder();
200 // quad-dominated surface meshing
204 mparams.quad = NETGENPlugin_Hypothesis_2D::GetDefaultQuadAllowed();
205 _fineness = NETGENPlugin_Hypothesis::GetDefaultFineness();
206 mparams.uselocalh = NETGENPlugin_Hypothesis::GetDefaultSurfaceCurvature();
207 netgen::merge_solids = NETGENPlugin_Hypothesis::GetDefaultFuseEdges();
210 //=============================================================================
214 //=============================================================================
216 void SetLocalSize(TopoDS_Shape GeomShape, double LocalSize)
218 if ( GeomShape.IsNull() ) return;
219 TopAbs_ShapeEnum GeomType = GeomShape.ShapeType();
220 if (GeomType == TopAbs_COMPOUND) {
221 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()) {
222 SetLocalSize(it.Value(), LocalSize);
227 if (! ShapesWithLocalSize.Contains(GeomShape))
228 key = ShapesWithLocalSize.Add(GeomShape);
230 key = ShapesWithLocalSize.FindIndex(GeomShape);
231 if (GeomType == TopAbs_VERTEX) {
232 VertexId2LocalSize[key] = LocalSize;
233 } else if (GeomType == TopAbs_EDGE) {
234 EdgeId2LocalSize[key] = LocalSize;
235 } else if (GeomType == TopAbs_FACE) {
236 FaceId2LocalSize[key] = LocalSize;
240 //=============================================================================
242 * Pass parameters to NETGEN
244 //=============================================================================
245 void NETGENPlugin_Mesher::SetParameters(const NETGENPlugin_Hypothesis* hyp)
249 netgen::MeshingParameters& mparams = netgen::mparam;
250 // Initialize global NETGEN parameters:
251 // maximal mesh segment size
252 mparams.maxh = hyp->GetMaxSize();
253 // maximal mesh element linear size
254 mparams.minh = hyp->GetMinSize();
255 // minimal number of segments per edge
256 mparams.segmentsperedge = hyp->GetNbSegPerEdge();
257 // rate of growth of size between elements
258 mparams.grading = hyp->GetGrowthRate();
259 // safety factor for curvatures (elements per radius)
260 mparams.curvaturesafety = hyp->GetNbSegPerRadius();
261 // create elements of second order
262 mparams.secondorder = hyp->GetSecondOrder() ? 1 : 0;
263 // quad-dominated surface meshing
264 // only triangles are allowed for volumic mesh (before realizing IMP 0021676)
266 mparams.quad = hyp->GetQuadAllowed() ? 1 : 0;
267 _optimize = hyp->GetOptimize();
268 _fineness = hyp->GetFineness();
269 mparams.uselocalh = hyp->GetSurfaceCurvature();
270 netgen::merge_solids = hyp->GetFuseEdges();
273 SMESH_Gen_i* smeshGen_i = SMESH_Gen_i::GetSMESHGen();
274 CORBA::Object_var anObject = smeshGen_i->GetNS()->Resolve("/myStudyManager");
275 SALOMEDS::StudyManager_var aStudyMgr = SALOMEDS::StudyManager::_narrow(anObject);
276 SALOMEDS::Study_var myStudy = aStudyMgr->GetStudyByID(hyp->GetStudyId());
278 const NETGENPlugin_Hypothesis::TLocalSize localSizes = hyp->GetLocalSizesAndEntries();
279 NETGENPlugin_Hypothesis::TLocalSize::const_iterator it = localSizes.begin();
280 for (it ; it != localSizes.end() ; it++)
282 std::string entry = (*it).first;
283 double val = (*it).second;
285 GEOM::GEOM_Object_var aGeomObj;
286 TopoDS_Shape S = TopoDS_Shape();
287 SALOMEDS::SObject_var aSObj = myStudy->FindObjectID( entry.c_str() );
288 if (!aSObj->_is_nil()) {
289 CORBA::Object_var obj = aSObj->GetObject();
290 aGeomObj = GEOM::GEOM_Object::_narrow(obj);
293 if ( !aGeomObj->_is_nil() )
294 S = smeshGen_i->GeomObjectToShape( aGeomObj.in() );
296 SetLocalSize(S, val);
301 //=============================================================================
303 * Pass simple parameters to NETGEN
305 //=============================================================================
307 void NETGENPlugin_Mesher::SetParameters(const NETGENPlugin_SimpleHypothesis_2D* hyp)
311 SetDefaultParameters();
314 //=============================================================================
316 * Link - a pair of integer numbers
318 //=============================================================================
322 Link(int _n1, int _n2) : n1(_n1), n2(_n2) {}
323 Link() : n1(0), n2(0) {}
326 int HashCode(const Link& aLink, int aLimit)
328 return HashCode(aLink.n1 + aLink.n2, aLimit);
331 Standard_Boolean IsEqual(const Link& aLink1, const Link& aLink2)
333 return (aLink1.n1 == aLink2.n1 && aLink1.n2 == aLink2.n2 ||
334 aLink1.n1 == aLink2.n2 && aLink1.n2 == aLink2.n1);
339 //================================================================================
341 * \brief return id of netgen point corresponding to SMDS node
343 //================================================================================
344 typedef map< const SMDS_MeshNode*, int > TNode2IdMap;
346 int ngNodeId( const SMDS_MeshNode* node,
347 netgen::Mesh& ngMesh,
348 TNode2IdMap& nodeNgIdMap)
350 int newNgId = ngMesh.GetNP() + 1;
352 TNode2IdMap::iterator node_id = nodeNgIdMap.insert( make_pair( node, newNgId )).first;
354 if ( node_id->second == newNgId)
356 #if defined(DUMP_SEGMENTS) || defined(DUMP_TRIANGLES)
357 cout << "Ng " << newNgId << " - " << node;
359 netgen::MeshPoint p( netgen::Point<3> (node->X(), node->Y(), node->Z()) );
360 ngMesh.AddPoint( p );
362 return node_id->second;
365 //================================================================================
367 * \brief Return computed EDGEs connected to the given one
369 //================================================================================
371 list< TopoDS_Edge > getConnectedEdges( const TopoDS_Edge& edge,
372 const TopoDS_Face& face,
373 const set< SMESH_subMesh* > & computedSM,
374 const SMESH_MesherHelper& helper,
375 map< SMESH_subMesh*, set< int > >& addedEdgeSM2Faces)
378 list< TopoDS_Edge > edges;
379 list< int > nbEdgesInWire;
380 int nbWires = SMESH_Block::GetOrderedEdges( face, edges, nbEdgesInWire);
382 // find <edge> within <edges>
383 list< TopoDS_Edge >::iterator eItFwd = edges.begin();
384 for ( ; eItFwd != edges.end(); ++eItFwd )
385 if ( edge.IsSame( *eItFwd ))
387 if ( eItFwd == edges.end()) return list< TopoDS_Edge>();
389 if ( eItFwd->Orientation() >= TopAbs_INTERNAL )
391 // connected INTERNAL edges returned from GetOrderedEdges() are wrongly oriented
392 // so treat each INTERNAL edge separately
393 TopoDS_Edge e = *eItFwd;
395 edges.push_back( e );
399 // get all computed EDGEs connected to <edge>
401 list< TopoDS_Edge >::iterator eItBack = eItFwd, ePrev;
402 TopoDS_Vertex vCommon;
403 TopTools_MapOfShape eAdded; // map used not to add a seam edge twice to <edges>
406 // put edges before <edge> to <edges> back
407 while ( edges.begin() != eItFwd )
408 edges.splice( edges.end(), edges, edges.begin() );
412 while ( ++eItFwd != edges.end() )
414 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( *eItFwd );
416 bool connected = TopExp::CommonVertex( *ePrev, *eItFwd, vCommon );
417 bool computed = sm->IsMeshComputed();
418 bool added = addedEdgeSM2Faces[ sm ].count( helper.GetSubShapeID() );
419 bool doubled = !eAdded.Add( *eItFwd );
420 bool orientOK = (( ePrev ->Orientation() < TopAbs_INTERNAL ) ==
421 ( eItFwd->Orientation() < TopAbs_INTERNAL ) );
422 if ( !connected || !computed || !orientOK || added || doubled )
424 // stop advancement; move edges from tail to head
425 while ( edges.back() != *ePrev )
426 edges.splice( edges.begin(), edges, --edges.end() );
432 while ( eItBack != edges.begin() )
436 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( *eItBack );
438 bool connected = TopExp::CommonVertex( *ePrev, *eItBack, vCommon );
439 bool computed = sm->IsMeshComputed();
440 bool added = addedEdgeSM2Faces[ sm ].count( helper.GetSubShapeID() );
441 bool doubled = !eAdded.Add( *eItBack );
442 bool orientOK = (( ePrev ->Orientation() < TopAbs_INTERNAL ) ==
443 ( eItBack->Orientation() < TopAbs_INTERNAL ) );
444 if ( !connected || !computed || !orientOK || added || doubled)
447 edges.erase( edges.begin(), ePrev );
451 if ( edges.front() != edges.back() )
453 // assure that the 1st vertex is meshed
454 TopoDS_Edge eLast = edges.back();
455 while ( !SMESH_Algo::VertexNode( SMESH_MesherHelper::IthVertex( 0, edges.front()), helper.GetMeshDS())
457 edges.front() != eLast )
458 edges.splice( edges.end(), edges, edges.begin() );
463 //================================================================================
465 * \brief Make triangulation of a shape precise enough
467 //================================================================================
469 void updateTriangulation( const TopoDS_Shape& shape )
471 // static set< Poly_Triangulation* > updated;
473 // TopLoc_Location loc;
474 // TopExp_Explorer fExp( shape, TopAbs_FACE );
475 // for ( ; fExp.More(); fExp.Next() )
477 // Handle(Poly_Triangulation) triangulation =
478 // BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
479 // if ( triangulation.IsNull() ||
480 // updated.insert( triangulation.operator->() ).second )
482 // BRepTools::Clean (shape);
485 BRepMesh_IncrementalMesh e(shape, 0.01, true);
487 catch (Standard_Failure)
490 // updated.erase( triangulation.operator->() );
491 // triangulation = BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
492 // updated.insert( triangulation.operator->() );
496 //================================================================================
498 * \brief Returns a medium node either existing in SMESH of created by NETGEN
499 * \param [in] corner1 - corner node 1
500 * \param [in] corner2 - corner node 2
501 * \param [in] defaultMedium - the node created by NETGEN
502 * \param [in] helper - holder of medium nodes existing in SMESH
503 * \return const SMDS_MeshNode* - the result node
505 //================================================================================
507 const SMDS_MeshNode* mediumNode( const SMDS_MeshNode* corner1,
508 const SMDS_MeshNode* corner2,
509 const SMDS_MeshNode* defaultMedium,
510 const SMESH_MesherHelper* helper)
514 TLinkNodeMap::const_iterator l2n =
515 helper->GetTLinkNodeMap().find( SMESH_TLink( corner1, corner2 ));
516 if ( l2n != helper->GetTLinkNodeMap().end() )
517 defaultMedium = l2n->second;
519 return defaultMedium;
522 //================================================================================
524 * \brief Assure that mesh on given shapes is quadratic
526 //================================================================================
528 void makeQuadratic( const TopTools_IndexedMapOfShape& shapes,
531 for ( int i = 1; i <= shapes.Extent(); ++i )
533 SMESHDS_SubMesh* smDS = mesh->GetMeshDS()->MeshElements( shapes(i) );
534 if ( !smDS ) continue;
535 SMDS_ElemIteratorPtr elemIt = smDS->GetElements();
536 if ( !elemIt->more() ) continue;
537 const SMDS_MeshElement* e = elemIt->next();
538 if ( !e || e->IsQuadratic() )
541 TIDSortedElemSet elems;
543 while ( elemIt->more() )
544 elems.insert( elems.end(), elemIt->next() );
546 SMESH_MeshEditor( mesh ).ConvertToQuadratic( /*3d=*/false, elems, /*biQuad=*/false );
552 //================================================================================
554 * \brief Initialize netgen::OCCGeometry with OCCT shape
556 //================================================================================
558 void NETGENPlugin_Mesher::PrepareOCCgeometry(netgen::OCCGeometry& occgeo,
559 const TopoDS_Shape& shape,
561 list< SMESH_subMesh* > * meshedSM,
562 NETGENPlugin_Internals* intern)
564 updateTriangulation( shape );
567 BRepBndLib::Add (shape, bb);
568 double x1,y1,z1,x2,y2,z2;
569 bb.Get (x1,y1,z1,x2,y2,z2);
570 MESSAGE("shape bounding box:\n" <<
571 "(" << x1 << " " << y1 << " " << z1 << ") " <<
572 "(" << x2 << " " << y2 << " " << z2 << ")");
573 netgen::Point<3> p1 = netgen::Point<3> (x1,y1,z1);
574 netgen::Point<3> p2 = netgen::Point<3> (x2,y2,z2);
575 occgeo.boundingbox = netgen::Box<3> (p1,p2);
577 occgeo.shape = shape;
580 // fill maps of shapes of occgeo with not yet meshed subshapes
582 // get root submeshes
583 list< SMESH_subMesh* > rootSM;
584 const int shapeID = mesh.GetMeshDS()->ShapeToIndex( shape );
585 if ( shapeID > 0 ) { // SMESH_subMesh with ID 0 may exist, don't use it!
586 rootSM.push_back( mesh.GetSubMesh( shape ));
589 for ( TopoDS_Iterator it( shape ); it.More(); it.Next() )
590 rootSM.push_back( mesh.GetSubMesh( it.Value() ));
593 // add subshapes of empty submeshes
594 list< SMESH_subMesh* >::iterator rootIt = rootSM.begin(), rootEnd = rootSM.end();
595 for ( ; rootIt != rootEnd; ++rootIt ) {
596 SMESH_subMesh * root = *rootIt;
597 SMESH_subMeshIteratorPtr smIt = root->getDependsOnIterator(/*includeSelf=*/true,
598 /*complexShapeFirst=*/true);
599 // to find a right orientation of subshapes (PAL20462)
600 TopTools_IndexedMapOfShape subShapes;
601 TopExp::MapShapes(root->GetSubShape(), subShapes);
602 while ( smIt->more() )
604 SMESH_subMesh* sm = smIt->next();
605 TopoDS_Shape shape = sm->GetSubShape();
606 if ( intern && intern->isShapeToPrecompute( shape ))
608 if ( !meshedSM || sm->IsEmpty() )
610 if ( shape.ShapeType() != TopAbs_VERTEX )
611 shape = subShapes( subShapes.FindIndex( shape ));// shape -> index -> oriented shape
612 if ( shape.Orientation() >= TopAbs_INTERNAL )
613 shape.Orientation( TopAbs_FORWARD ); // isuue 0020676
614 switch ( shape.ShapeType() ) {
615 case TopAbs_FACE : occgeo.fmap.Add( shape ); break;
616 case TopAbs_EDGE : occgeo.emap.Add( shape ); break;
617 case TopAbs_VERTEX: occgeo.vmap.Add( shape ); break;
618 case TopAbs_SOLID :occgeo.somap.Add( shape ); break;
622 // collect submeshes of meshed shapes
625 const int dim = SMESH_Gen::GetShapeDim( shape );
626 meshedSM[ dim ].push_back( sm );
630 occgeo.facemeshstatus.SetSize (occgeo.fmap.Extent());
631 occgeo.facemeshstatus = 0;
632 occgeo.face_maxh_modified.SetSize(occgeo.fmap.Extent());
633 occgeo.face_maxh_modified = 0;
634 occgeo.face_maxh.SetSize(occgeo.fmap.Extent());
635 occgeo.face_maxh = netgen::mparam.maxh;
638 //================================================================================
640 * \brief Return a default min size value suitable for the given geometry.
642 //================================================================================
644 double NETGENPlugin_Mesher::GetDefaultMinSize(const TopoDS_Shape& geom,
645 const double maxSize)
647 updateTriangulation( geom );
651 const int* pi[4] = { &i1, &i2, &i3, &i1 };
654 TopExp_Explorer fExp( geom, TopAbs_FACE );
655 for ( ; fExp.More(); fExp.Next() )
657 Handle(Poly_Triangulation) triangulation =
658 BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
659 if ( triangulation.IsNull() ) continue;
660 const double fTol = BRep_Tool::Tolerance( TopoDS::Face( fExp.Current() ));
661 const TColgp_Array1OfPnt& points = triangulation->Nodes();
662 const Poly_Array1OfTriangle& trias = triangulation->Triangles();
663 for ( int iT = trias.Lower(); iT <= trias.Upper(); ++iT )
665 trias(iT).Get( i1, i2, i3 );
666 for ( int j = 0; j < 3; ++j )
668 double dist2 = points(*pi[j]).SquareDistance( points( *pi[j+1] ));
669 if ( dist2 < minh && fTol*fTol < dist2 )
671 bb.Add( points(*pi[j]));
675 if ( minh > 0.25 * bb.SquareExtent() ) // simple geometry, rough triangulation
677 minh = 1e-3 * sqrt( bb.SquareExtent());
678 //cout << "BND BOX minh = " <<minh << endl;
682 minh = 3 * sqrt( minh ); // triangulation for visualization is rather fine
683 //cout << "TRIANGULATION minh = " <<minh << endl;
685 if ( minh > 0.5 * maxSize )
691 //================================================================================
693 * \brief Restrict size of elements at a given point
695 //================================================================================
697 void NETGENPlugin_Mesher::RestrictLocalSize(netgen::Mesh& ngMesh,
700 const bool overrideMinH)
702 if ( overrideMinH && netgen::mparam.minh > size )
704 ngMesh.SetMinimalH( size );
705 netgen::mparam.minh = size;
707 netgen::Point3d pi(p.X(), p.Y(), p.Z());
708 ngMesh.RestrictLocalH( pi, size );
711 //================================================================================
713 * \brief fill ngMesh with nodes and elements of computed submeshes
715 //================================================================================
717 bool NETGENPlugin_Mesher::FillNgMesh(netgen::OCCGeometry& occgeom,
718 netgen::Mesh& ngMesh,
719 vector<const SMDS_MeshNode*>& nodeVec,
720 const list< SMESH_subMesh* > & meshedSM,
721 SMESH_MesherHelper* quadHelper,
722 SMESH_ProxyMesh::Ptr proxyMesh)
724 TNode2IdMap nodeNgIdMap;
725 for ( int i = 1; i < nodeVec.size(); ++i )
726 nodeNgIdMap.insert( make_pair( nodeVec[i], i ));
728 TopTools_MapOfShape visitedShapes;
729 map< SMESH_subMesh*, set< int > > visitedEdgeSM2Faces;
730 set< SMESH_subMesh* > computedSM( meshedSM.begin(), meshedSM.end() );
732 SMESH_MesherHelper helper (*_mesh);
734 int faceNgID = ngMesh.GetNFD();
736 list< SMESH_subMesh* >::const_iterator smIt, smEnd = meshedSM.end();
737 for ( smIt = meshedSM.begin(); smIt != smEnd; ++smIt )
739 SMESH_subMesh* sm = *smIt;
740 if ( !visitedShapes.Add( sm->GetSubShape() ))
743 const SMESHDS_SubMesh * smDS = sm->GetSubMeshDS();
744 if ( !smDS ) continue;
746 switch ( sm->GetSubShape().ShapeType() )
748 case TopAbs_EDGE: { // EDGE
749 // ----------------------
750 TopoDS_Edge geomEdge = TopoDS::Edge( sm->GetSubShape() );
751 if ( geomEdge.Orientation() >= TopAbs_INTERNAL )
752 geomEdge.Orientation( TopAbs_FORWARD ); // issue 0020676
754 // Add ng segments for each not meshed FACE the EDGE bounds
755 PShapeIteratorPtr fIt = helper.GetAncestors( geomEdge, *sm->GetFather(), TopAbs_FACE );
756 while ( const TopoDS_Shape * anc = fIt->next() )
758 faceNgID = occgeom.fmap.FindIndex( *anc );
760 continue; // meshed face
762 int faceSMDSId = helper.GetMeshDS()->ShapeToIndex( *anc );
763 if ( visitedEdgeSM2Faces[ sm ].count( faceSMDSId ))
764 continue; // already treated EDGE
766 TopoDS_Face face = TopoDS::Face( occgeom.fmap( faceNgID ));
767 if ( face.Orientation() >= TopAbs_INTERNAL )
768 face.Orientation( TopAbs_FORWARD ); // issue 0020676
770 // get all meshed EDGEs of the FACE connected to geomEdge (issue 0021140)
771 helper.SetSubShape( face );
772 list< TopoDS_Edge > edges = getConnectedEdges( geomEdge, face, computedSM, helper,
773 visitedEdgeSM2Faces );
775 continue; // wrong ancestor?
777 // find out orientation of <edges> within <face>
778 TopoDS_Edge eNotSeam = edges.front();
779 if ( helper.HasSeam() )
781 list< TopoDS_Edge >::iterator eIt = edges.begin();
782 while ( helper.IsRealSeam( *eIt )) ++eIt;
783 if ( eIt != edges.end() )
786 TopAbs_Orientation fOri = helper.GetSubShapeOri( face, eNotSeam );
787 bool isForwad = ( fOri == eNotSeam.Orientation() || fOri >= TopAbs_INTERNAL );
789 // get all nodes from connected <edges>
790 bool isQuad = smDS->NbElements() ? smDS->GetElements()->next()->IsQuadratic() : false;
791 StdMeshers_FaceSide fSide( face, edges, _mesh, isForwad, isQuad );
792 const vector<UVPtStruct>& points = fSide.GetUVPtStruct();
793 int i, nbSeg = fSide.NbSegments();
795 // remember EDGEs of fSide to treat only once
796 for ( int iE = 0; iE < fSide.NbEdges(); ++iE )
797 visitedEdgeSM2Faces[ helper.GetMesh()->GetSubMesh( fSide.Edge(iE )) ].insert(faceSMDSId);
799 double otherSeamParam = 0;
804 int prevNgId = ngNodeId( points[0].node, ngMesh, nodeNgIdMap );
806 for ( i = 0; i < nbSeg; ++i )
808 const UVPtStruct& p1 = points[ i ];
809 const UVPtStruct& p2 = points[ i+1 ];
811 if ( p1.node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX ) //an EDGE begins
814 if ( helper.IsRealSeam( p1.node->getshapeId() ))
816 TopoDS_Edge e = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
817 isSeam = helper.IsRealSeam( e );
820 otherSeamParam = helper.GetOtherParam( helper.GetPeriodicIndex() & 1 ? p2.u : p2.v );
827 seg[1] = prevNgId = ngNodeId( p2.node, ngMesh, nodeNgIdMap );
828 // node param on curve
829 seg.epgeominfo[ 0 ].dist = p1.param;
830 seg.epgeominfo[ 1 ].dist = p2.param;
832 seg.epgeominfo[ 0 ].u = p1.u;
833 seg.epgeominfo[ 0 ].v = p1.v;
834 seg.epgeominfo[ 1 ].u = p2.u;
835 seg.epgeominfo[ 1 ].v = p2.v;
837 //geomEdge = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
838 //seg.epgeominfo[ 0 ].edgenr = seg.epgeominfo[ 1 ].edgenr = occgeom.emap.FindIndex( geomEdge );
840 //seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
841 seg.si = faceNgID; // = geom.fmap.FindIndex (face);
842 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
843 ngMesh.AddSegment (seg);
845 SMESH_TNodeXYZ np1( p1.node ), np2( p2.node );
846 RestrictLocalSize( ngMesh, 0.5*(np1+np2), (np1-np2).Modulus() );
849 cout << "Segment: " << seg.edgenr << " on SMESH face " << helper.GetMeshDS()->ShapeToIndex( face ) << endl
850 << "\tface index: " << seg.si << endl
851 << "\tp1: " << seg[0] << endl
852 << "\tp2: " << seg[1] << endl
853 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
854 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
855 //<< "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
856 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
857 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl;
858 //<< "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
862 if ( helper.GetPeriodicIndex() && 1 ) {
863 seg.epgeominfo[ 0 ].u = otherSeamParam;
864 seg.epgeominfo[ 1 ].u = otherSeamParam;
865 swap (seg.epgeominfo[0].v, seg.epgeominfo[1].v);
867 seg.epgeominfo[ 0 ].v = otherSeamParam;
868 seg.epgeominfo[ 1 ].v = otherSeamParam;
869 swap (seg.epgeominfo[0].u, seg.epgeominfo[1].u);
871 swap (seg[0], seg[1]);
872 swap (seg.epgeominfo[0].dist, seg.epgeominfo[1].dist);
873 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
874 ngMesh.AddSegment (seg);
876 cout << "Segment: " << seg.edgenr << endl
877 << "\t is SEAM (reverse) of the previous. "
878 << " Other " << (helper.GetPeriodicIndex() && 1 ? "U" : "V")
879 << " = " << otherSeamParam << endl;
882 else if ( fOri == TopAbs_INTERNAL )
884 swap (seg[0], seg[1]);
885 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
886 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
887 ngMesh.AddSegment (seg);
889 cout << "Segment: " << seg.edgenr << endl << "\t is REVERSE of the previous" << endl;
893 } // loop on geomEdge ancestors
895 if ( quadHelper ) // remember medium nodes of sub-meshes
897 SMDS_ElemIteratorPtr edges = smDS->GetElements();
898 while ( edges->more() )
900 const SMDS_MeshElement* e = edges->next();
901 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshEdge*>( e )))
907 } // case TopAbs_EDGE
909 case TopAbs_FACE: { // FACE
910 // ----------------------
911 const TopoDS_Face& geomFace = TopoDS::Face( sm->GetSubShape() );
912 helper.SetSubShape( geomFace );
913 bool isInternalFace = ( geomFace.Orientation() == TopAbs_INTERNAL );
915 // Find solids the geomFace bounds
916 int solidID1 = 0, solidID2 = 0;
917 StdMeshers_QuadToTriaAdaptor* quadAdaptor =
918 dynamic_cast<StdMeshers_QuadToTriaAdaptor*>( proxyMesh.get() );
921 solidID1 = occgeom.somap.FindIndex( quadAdaptor->GetShape() );
925 PShapeIteratorPtr solidIt = helper.GetAncestors( geomFace, *sm->GetFather(), TopAbs_SOLID);
926 while ( const TopoDS_Shape * solid = solidIt->next() )
928 int id = occgeom.somap.FindIndex ( *solid );
929 if ( solidID1 && id != solidID1 ) solidID2 = id;
933 // Add ng face descriptors of meshed faces
935 ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(faceNgID, solidID1, solidID2, 0));
937 // if second oreder is required, even already meshed faces must be passed to NETGEN
938 int fID = occgeom.fmap.Add( geomFace );
939 while ( fID < faceNgID ) // geomFace is already in occgeom.fmap, add a copy
940 fID = occgeom.fmap.Add( BRepBuilderAPI_Copy( geomFace, /*copyGeom=*/false ));
941 // Problem with the second order in a quadrangular mesh remains.
942 // 1) All quadrangles generated by NETGEN are moved to an inexistent face
943 // by FillSMesh() (find "AddFaceDescriptor")
944 // 2) Temporary triangles generated by StdMeshers_QuadToTriaAdaptor
945 // are on faces where quadrangles were.
946 // Due to these 2 points, wrong geom faces are used while conversion to qudratic
947 // of the mentioned above quadrangles and triangles
949 // Orient the face correctly in solidID1 (issue 0020206)
950 bool reverse = false;
952 TopoDS_Shape solid = occgeom.somap( solidID1 );
953 TopAbs_Orientation faceOriInSolid = helper.GetSubShapeOri( solid, geomFace );
954 if ( faceOriInSolid >= 0 )
956 helper.IsReversedSubMesh( TopoDS::Face( geomFace.Oriented( faceOriInSolid )));
959 // Add surface elements
961 netgen::Element2d tri(3);
962 tri.SetIndex ( faceNgID );
965 #ifdef DUMP_TRIANGLES
966 cout << "SMESH face " << helper.GetMeshDS()->ShapeToIndex( geomFace )
967 << " internal="<<isInternalFace << endl;
970 smDS = proxyMesh->GetSubMesh( geomFace );
972 SMDS_ElemIteratorPtr faces = smDS->GetElements();
973 while ( faces->more() )
975 const SMDS_MeshElement* f = faces->next();
976 if ( f->NbNodes() % 3 != 0 ) // not triangle
978 PShapeIteratorPtr solidIt=helper.GetAncestors(geomFace,*sm->GetFather(),TopAbs_SOLID);
979 if ( const TopoDS_Shape * solid = solidIt->next() )
980 sm = _mesh->GetSubMesh( *solid );
981 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
982 smError.reset( new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH,"Not triangle submesh"));
983 smError->myBadElements.push_back( f );
987 for ( int i = 0; i < 3; ++i )
989 const SMDS_MeshNode* node = f->GetNode( i ), * inFaceNode=0;
991 // get node UV on face
992 int shapeID = node->getshapeId();
993 if ( helper.IsSeamShape( shapeID ))
994 if ( helper.IsSeamShape( f->GetNodeWrap( i+1 )->getshapeId() ))
995 inFaceNode = f->GetNodeWrap( i-1 );
997 inFaceNode = f->GetNodeWrap( i+1 );
998 gp_XY uv = helper.GetNodeUV( geomFace, node, inFaceNode );
1000 int ind = reverse ? 3-i : i+1;
1001 tri.GeomInfoPi(ind).u = uv.X();
1002 tri.GeomInfoPi(ind).v = uv.Y();
1003 tri.PNum (ind) = ngNodeId( node, ngMesh, nodeNgIdMap );
1006 ngMesh.AddSurfaceElement (tri);
1007 #ifdef DUMP_TRIANGLES
1008 cout << tri << endl;
1011 if ( isInternalFace )
1013 swap( tri[1], tri[2] );
1014 ngMesh.AddSurfaceElement (tri);
1015 #ifdef DUMP_TRIANGLES
1016 cout << tri << endl;
1021 if ( quadHelper ) // remember medium nodes of sub-meshes
1023 SMDS_ElemIteratorPtr faces = smDS->GetElements();
1024 while ( faces->more() )
1026 const SMDS_MeshElement* f = faces->next();
1027 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshFace*>( f )))
1033 } // case TopAbs_FACE
1035 case TopAbs_VERTEX: { // VERTEX
1036 // --------------------------
1037 // issue 0021405. Add node only if a VERTEX is shared by a not meshed EDGE,
1038 // else netgen removes a free node and nodeVector becomes invalid
1039 PShapeIteratorPtr ansIt = helper.GetAncestors( sm->GetSubShape(),
1043 while ( const TopoDS_Shape* e = ansIt->next() )
1045 SMESH_subMesh* eSub = helper.GetMesh()->GetSubMesh( *e );
1046 if (( toAdd = ( eSub->IsEmpty() && !SMESH_Algo::isDegenerated( TopoDS::Edge( *e )))))
1051 SMDS_NodeIteratorPtr nodeIt = smDS->GetNodes();
1052 if ( nodeIt->more() )
1053 ngNodeId( nodeIt->next(), ngMesh, nodeNgIdMap );
1059 } // loop on submeshes
1062 nodeVec.resize( ngMesh.GetNP() + 1 );
1063 TNode2IdMap::iterator node_NgId, nodeNgIdEnd = nodeNgIdMap.end();
1064 for ( node_NgId = nodeNgIdMap.begin(); node_NgId != nodeNgIdEnd; ++node_NgId)
1065 nodeVec[ node_NgId->second ] = node_NgId->first;
1070 //================================================================================
1072 * \brief Duplicate mesh faces on internal geom faces
1074 //================================================================================
1076 void NETGENPlugin_Mesher::FixIntFaces(const netgen::OCCGeometry& occgeom,
1077 netgen::Mesh& ngMesh,
1078 NETGENPlugin_Internals& internalShapes)
1080 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1082 // find ng indices of internal faces
1084 for ( int ngFaceID = 1; ngFaceID <= occgeom.fmap.Extent(); ++ngFaceID )
1086 int smeshID = meshDS->ShapeToIndex( occgeom.fmap( ngFaceID ));
1087 if ( internalShapes.isInternalShape( smeshID ))
1088 ngFaceIds.insert( ngFaceID );
1090 if ( !ngFaceIds.empty() )
1093 int i, nbFaces = ngMesh.GetNSE();
1094 for (int i = 1; i <= nbFaces; ++i)
1096 netgen::Element2d elem = ngMesh.SurfaceElement(i);
1097 if ( ngFaceIds.count( elem.GetIndex() ))
1099 swap( elem[1], elem[2] );
1100 ngMesh.AddSurfaceElement (elem);
1108 //================================================================================
1109 // define gp_XY_Subtracted pointer to function calling gp_XY::Subtracted(gp_XY)
1110 gp_XY_FunPtr(Subtracted);
1111 //gp_XY_FunPtr(Added);
1113 //================================================================================
1115 * \brief Evaluate distance between two 2d points along the surface
1117 //================================================================================
1119 double evalDist( const gp_XY& uv1,
1121 const Handle(Geom_Surface)& surf,
1122 const int stopHandler=-1)
1124 if ( stopHandler > 0 ) // continue recursion
1126 gp_XY mid = SMESH_MesherHelper::GetMiddleUV( surf, uv1, uv2 );
1127 return evalDist( uv1,mid, surf, stopHandler-1 ) + evalDist( mid,uv2, surf, stopHandler-1 );
1129 double dist3D = surf->Value( uv1.X(), uv1.Y() ).Distance( surf->Value( uv2.X(), uv2.Y() ));
1130 if ( stopHandler == 0 ) // stop recursion
1133 // start recursion if necessary
1134 double dist2D = SMESH_MesherHelper::applyIn2D(surf, uv1, uv2, gp_XY_Subtracted, 0).Modulus();
1135 if ( fabs( dist3D - dist2D ) < dist2D * 1e-10 )
1136 return dist3D; // equal parametrization of a planar surface
1138 return evalDist( uv1, uv2, surf, 3 ); // start recursion
1141 //================================================================================
1143 * \brief Data of vertex internal in geom face
1145 //================================================================================
1149 gp_XY uv; //!< UV in face parametric space
1150 int ngId; //!< ng id of corrsponding node
1151 gp_XY uvClose; //!< UV of closest boundary node
1152 int ngIdClose; //!< ng id of closest boundary node
1155 //================================================================================
1157 * \brief Data of vertex internal in solid
1159 //================================================================================
1163 int ngId; //!< ng id of corresponding node
1164 int ngIdClose; //!< ng id of closest 2d mesh element
1165 int ngIdCloseN; //!< ng id of closest node of the closest 2d mesh element
1168 inline double dist2(const netgen::MeshPoint& p1, const netgen::MeshPoint& p2)
1170 return gp_Pnt( NGPOINT_COORDS(p1)).SquareDistance( gp_Pnt( NGPOINT_COORDS(p2)));
1174 //================================================================================
1176 * \brief Make netgen take internal vertices in faces into account by adding
1177 * segments including internal vertices
1179 * This function works in supposition that 1D mesh is already computed in ngMesh
1181 //================================================================================
1183 void NETGENPlugin_Mesher::AddIntVerticesInFaces(const netgen::OCCGeometry& occgeom,
1184 netgen::Mesh& ngMesh,
1185 vector<const SMDS_MeshNode*>& nodeVec,
1186 NETGENPlugin_Internals& internalShapes)
1188 if ( nodeVec.size() < ngMesh.GetNP() )
1189 nodeVec.resize( ngMesh.GetNP(), 0 );
1191 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1192 SMESH_MesherHelper helper( internalShapes.getMesh() );
1194 const map<int,list<int> >& face2Vert = internalShapes.getFacesWithVertices();
1195 map<int,list<int> >::const_iterator f2v = face2Vert.begin();
1196 for ( ; f2v != face2Vert.end(); ++f2v )
1198 const TopoDS_Face& face = TopoDS::Face( meshDS->IndexToShape( f2v->first ));
1199 if ( face.IsNull() ) continue;
1200 int faceNgID = occgeom.fmap.FindIndex (face);
1201 if ( faceNgID < 0 ) continue;
1203 TopLoc_Location loc;
1204 Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc);
1206 helper.SetSubShape( face );
1207 helper.SetElementsOnShape( true );
1209 // Get data of internal vertices and add them to ngMesh
1211 multimap< double, TIntVData > dist2VData; // sort vertices by distance from boundary nodes
1213 int i, nbSegInit = ngMesh.GetNSeg();
1215 // boundary characteristics
1216 double totSegLen2D = 0;
1219 const list<int>& iVertices = f2v->second;
1220 list<int>::const_iterator iv = iVertices.begin();
1221 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1224 // get node on vertex
1225 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1226 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1229 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1230 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1231 nV = SMESH_Algo::VertexNode( V, meshDS );
1232 if ( !nV ) continue;
1235 netgen::MeshPoint mp( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1236 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1237 vData.ngId = ngMesh.GetNP();
1238 nodeVec.push_back( nV );
1242 vData.uv = helper.GetNodeUV( face, nV, 0, &uvOK );
1243 if ( !uvOK ) helper.CheckNodeUV( face, nV, vData.uv, BRep_Tool::Tolerance(V),/*force=*/1);
1245 // loop on all segments of the face to find the node closest to vertex and to count
1246 // average segment 2d length
1247 double closeDist2 = numeric_limits<double>::max(), dist2;
1249 for (i = 1; i <= ngMesh.GetNSeg(); ++i)
1251 netgen::Segment & seg = ngMesh.LineSegment(i);
1252 if ( seg.si != faceNgID ) continue;
1254 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1256 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1257 if ( ngIdLast == seg[ iEnd ] ) continue;
1258 dist2 = helper.applyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1259 if ( dist2 < closeDist2 )
1260 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1261 ngIdLast = seg[ iEnd ];
1265 totSegLen2D += helper.applyIn2D(surf, uv[0], uv[1], gp_XY_Subtracted, false).Modulus();
1269 dist2VData.insert( make_pair( closeDist2, vData ));
1272 if ( totNbSeg == 0 ) break;
1273 double avgSegLen2d = totSegLen2D / totNbSeg;
1275 // Loop on vertices to add segments
1277 multimap< double, TIntVData >::iterator dist_vData = dist2VData.begin();
1278 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1280 double closeDist2 = dist_vData->first, dist2;
1281 TIntVData & vData = dist_vData->second;
1283 // try to find more close node among segments added for internal vertices
1284 for (i = nbSegInit+1; i <= ngMesh.GetNSeg(); ++i)
1286 netgen::Segment & seg = ngMesh.LineSegment(i);
1287 if ( seg.si != faceNgID ) continue;
1289 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1291 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1292 dist2 = helper.applyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1293 if ( dist2 < closeDist2 )
1294 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1297 // decide whether to use the closest node as the second end of segment or to
1298 // create a new point
1299 int segEnd1 = vData.ngId;
1300 int segEnd2 = vData.ngIdClose; // to use closest node
1301 gp_XY uvV = vData.uv, uvP = vData.uvClose;
1302 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1303 double nodeDist2D = sqrt( closeDist2 );
1304 double nodeDist3D = evalDist( vData.uv, vData.uvClose, surf );
1305 bool avgLenOK = ( avgSegLen2d < 0.75 * nodeDist2D );
1306 bool hintLenOK = ( segLenHint < 0.75 * nodeDist3D );
1307 //cout << "uvV " << uvV.X() <<","<<uvV.Y() << " ";
1308 if ( hintLenOK || avgLenOK )
1310 // create a point between the closest node and V
1313 double r = min( 0.5, ( hintLenOK ? segLenHint/nodeDist3D : avgSegLen2d/nodeDist2D ));
1314 // direction from V to closet node in 2D
1315 gp_Dir2d v2n( helper.applyIn2D(surf, uvP, uvV, gp_XY_Subtracted, false ));
1317 uvP = vData.uv + r * nodeDist2D * v2n.XY();
1318 gp_Pnt P = surf->Value( uvP.X(), uvP.Y() ).Transformed( loc );
1320 netgen::MeshPoint mp( netgen::Point<3> (P.X(), P.Y(), P.Z()));
1321 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1322 segEnd2 = ngMesh.GetNP();
1323 //cout << "Middle " << r << " uv " << uvP.X() << "," << uvP.Y() << "( " << ngMesh.Point(segEnd2).X()<<","<<ngMesh.Point(segEnd2).Y()<<","<<ngMesh.Point(segEnd2).Z()<<" )"<< endl;
1324 SMDS_MeshNode * nP = helper.AddNode(P.X(), P.Y(), P.Z());
1325 nodeVec.push_back( nP );
1327 //else cout << "at Node " << " uv " << uvP.X() << "," << uvP.Y() << endl;
1330 netgen::Segment seg;
1332 if ( segEnd1 > segEnd2 ) swap( segEnd1, segEnd2 ), swap( uvV, uvP );
1333 seg[0] = segEnd1; // ng node id
1334 seg[1] = segEnd2; // ng node id
1335 seg.edgenr = ngMesh.GetNSeg() + 1;// segment id
1338 seg.epgeominfo[ 0 ].dist = 0; // param on curve
1339 seg.epgeominfo[ 0 ].u = uvV.X();
1340 seg.epgeominfo[ 0 ].v = uvV.Y();
1341 seg.epgeominfo[ 1 ].dist = 1; // param on curve
1342 seg.epgeominfo[ 1 ].u = uvP.X();
1343 seg.epgeominfo[ 1 ].v = uvP.Y();
1345 // seg.epgeominfo[ 0 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1346 // seg.epgeominfo[ 1 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1348 ngMesh.AddSegment (seg);
1350 // add reverse segment
1351 swap (seg[0], seg[1]);
1352 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1353 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1354 ngMesh.AddSegment (seg);
1360 //================================================================================
1362 * \brief Make netgen take internal vertices in solids into account by adding
1363 * faces including internal vertices
1365 * This function works in supposition that 2D mesh is already computed in ngMesh
1367 //================================================================================
1369 void NETGENPlugin_Mesher::AddIntVerticesInSolids(const netgen::OCCGeometry& occgeom,
1370 netgen::Mesh& ngMesh,
1371 vector<const SMDS_MeshNode*>& nodeVec,
1372 NETGENPlugin_Internals& internalShapes)
1374 #ifdef DUMP_TRIANGLES_SCRIPT
1375 // create a python script making a mesh containing triangles added for internal vertices
1376 ofstream py(DUMP_TRIANGLES_SCRIPT);
1377 py << "import SMESH"<< endl
1378 << "from salome.smesh import smeshBuilder"<<endl
1379 << "smesh = smeshBuilder.New(salome.myStudy)"<<endl
1380 << "m = smesh.Mesh(name='triangles')" << endl;
1382 if ( nodeVec.size() < ngMesh.GetNP() )
1383 nodeVec.resize( ngMesh.GetNP(), 0 );
1385 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1386 SMESH_MesherHelper helper( internalShapes.getMesh() );
1388 const map<int,list<int> >& so2Vert = internalShapes.getSolidsWithVertices();
1389 map<int,list<int> >::const_iterator s2v = so2Vert.begin();
1390 for ( ; s2v != so2Vert.end(); ++s2v )
1392 const TopoDS_Shape& solid = meshDS->IndexToShape( s2v->first );
1393 if ( solid.IsNull() ) continue;
1394 int solidNgID = occgeom.somap.FindIndex (solid);
1395 if ( solidNgID < 0 && !occgeom.somap.IsEmpty() ) continue;
1397 helper.SetSubShape( solid );
1398 helper.SetElementsOnShape( true );
1400 // find ng indices of faces within the solid
1402 for (TopExp_Explorer fExp(solid, TopAbs_FACE); fExp.More(); fExp.Next() )
1403 ngFaceIds.insert( occgeom.fmap.FindIndex( fExp.Current() ));
1404 if ( ngFaceIds.size() == 1 && *ngFaceIds.begin() == 0 )
1405 ngFaceIds.insert( 1 );
1407 // Get data of internal vertices and add them to ngMesh
1409 multimap< double, TIntVSoData > dist2VData; // sort vertices by distance from ng faces
1411 int i, nbFaceInit = ngMesh.GetNSE();
1413 // boundary characteristics
1414 double totSegLen = 0;
1417 const list<int>& iVertices = s2v->second;
1418 list<int>::const_iterator iv = iVertices.begin();
1419 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1422 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1424 // get node on vertex
1425 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1428 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1429 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1430 nV = SMESH_Algo::VertexNode( V, meshDS );
1431 if ( !nV ) continue;
1434 netgen::MeshPoint mpV( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1435 ngMesh.AddPoint ( mpV, 1, netgen::FIXEDPOINT );
1436 vData.ngId = ngMesh.GetNP();
1437 nodeVec.push_back( nV );
1439 // loop on all 2d elements to find the one closest to vertex and to count
1440 // average segment length
1441 double closeDist2 = numeric_limits<double>::max(), avgDist2;
1442 for (i = 1; i <= ngMesh.GetNSE(); ++i)
1444 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1445 if ( !ngFaceIds.count( elem.GetIndex() )) continue;
1447 multimap< double, int> dist2nID; // sort nodes of element by distance from V
1448 for ( int j = 0; j < elem.GetNP(); ++j)
1450 netgen::MeshPoint mp = ngMesh.Point( elem[j] );
1451 double d2 = dist2( mpV, mp );
1452 dist2nID.insert( make_pair( d2, elem[j] ));
1453 avgDist2 += d2 / elem.GetNP();
1455 totNbSeg++, totSegLen+= sqrt( dist2( mp, ngMesh.Point( elem[(j+1)%elem.GetNP()])));
1457 double dist = dist2nID.begin()->first; //avgDist2;
1458 if ( dist < closeDist2 )
1459 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= dist;
1461 dist2VData.insert( make_pair( closeDist2, vData ));
1464 if ( totNbSeg == 0 ) break;
1465 double avgSegLen = totSegLen / totNbSeg;
1467 // Loop on vertices to add triangles
1469 multimap< double, TIntVSoData >::iterator dist_vData = dist2VData.begin();
1470 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1472 double closeDist2 = dist_vData->first;
1473 TIntVSoData & vData = dist_vData->second;
1475 const netgen::MeshPoint& mpV = ngMesh.Point( vData.ngId );
1477 // try to find more close face among ones added for internal vertices
1478 for (i = nbFaceInit+1; i <= ngMesh.GetNSE(); ++i)
1480 double avgDist2 = 0;
1481 multimap< double, int> dist2nID;
1482 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1483 for ( int j = 0; j < elem.GetNP(); ++j)
1485 double d = dist2( mpV, ngMesh.Point( elem[j] ));
1486 dist2nID.insert( make_pair( d, elem[j] ));
1487 avgDist2 += d / elem.GetNP();
1488 if ( avgDist2 < closeDist2 )
1489 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= avgDist2;
1492 // sort nodes of the closest face by angle with vector from V to the closest node
1493 const double tol = numeric_limits<double>::min();
1494 map< double, int > angle2ID;
1495 const netgen::Element2d& closeFace = ngMesh.SurfaceElement( vData.ngIdClose );
1496 netgen::MeshPoint mp[2];
1497 mp[0] = ngMesh.Point( vData.ngIdCloseN );
1498 gp_XYZ p1( NGPOINT_COORDS( mp[0] ));
1499 gp_XYZ pV( NGPOINT_COORDS( mpV ));
1500 gp_Vec v2p1( pV, p1 );
1501 double distN1 = v2p1.Magnitude();
1502 if ( distN1 <= tol ) continue;
1504 for ( int j = 0; j < closeFace.GetNP(); ++j)
1506 mp[1] = ngMesh.Point( closeFace[j] );
1507 gp_Vec v2p( pV, gp_Pnt( NGPOINT_COORDS( mp[1] )) );
1508 angle2ID.insert( make_pair( v2p1.Angle( v2p ), closeFace[j]));
1510 // get node with angle of 60 degrees or greater
1511 map< double, int >::iterator angle_id = angle2ID.lower_bound( 60. * M_PI / 180. );
1512 if ( angle_id == angle2ID.end() ) angle_id = --angle2ID.end();
1513 const double minAngle = 30. * M_PI / 180.;
1514 const double angle = angle_id->first;
1515 bool angleOK = ( angle > minAngle );
1517 // find points to create a triangle
1518 netgen::Element2d tri(3);
1520 tri[0] = vData.ngId;
1521 tri[1] = vData.ngIdCloseN; // to use the closest nodes
1522 tri[2] = angle_id->second; // to use the node with best angle
1524 // decide whether to use the closest node and the node with best angle or to create new ones
1525 for ( int isBestAngleN = 0; isBestAngleN < 2; ++isBestAngleN )
1527 bool createNew = !angleOK, distOK = true;
1529 int triInd = isBestAngleN ? 2 : 1;
1530 mp[isBestAngleN] = ngMesh.Point( tri[triInd] );
1535 double distN2 = sqrt( dist2( mpV, mp[isBestAngleN]));
1536 createNew = ( fabs( distN2 - distN1 ) > 0.25 * distN1 );
1538 else if ( angle < tol )
1540 v2p1.SetX( v2p1.X() + 1e-3 );
1546 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1547 bool avgLenOK = ( avgSegLen < 0.75 * distN1 );
1548 bool hintLenOK = ( segLenHint < 0.75 * distN1 );
1549 createNew = (createNew || avgLenOK || hintLenOK );
1550 // we create a new node not closer than 0.5 to the closest face
1551 // in order not to clash with other close face
1552 double r = min( 0.5, ( hintLenOK ? segLenHint : avgSegLen ) / distN1 );
1553 distFromV = r * distN1;
1557 // create a new point, between the node and the vertex if angleOK
1558 gp_XYZ p( NGPOINT_COORDS( mp[isBestAngleN] ));
1559 gp_Vec v2p( pV, p ); v2p.Normalize();
1560 if ( isBestAngleN && !angleOK )
1561 p = p1 + gp_Dir( v2p.XYZ() - v2p1.XYZ()).XYZ() * distN1 * 0.95;
1563 p = pV + v2p.XYZ() * distFromV;
1565 if ( !isBestAngleN ) p1 = p, distN1 = distFromV;
1567 mp[isBestAngleN].SetPoint( netgen::Point<3> (p.X(), p.Y(), p.Z()));
1568 ngMesh.AddPoint ( mp[isBestAngleN], 1, netgen::SURFACEPOINT );
1569 tri[triInd] = ngMesh.GetNP();
1570 nodeVec.push_back( helper.AddNode( p.X(), p.Y(), p.Z()) );
1573 ngMesh.AddSurfaceElement (tri);
1574 swap( tri[1], tri[2] );
1575 ngMesh.AddSurfaceElement (tri);
1577 #ifdef DUMP_TRIANGLES_SCRIPT
1578 py << "n1 = m.AddNode( "<< mpV(0)<<", "<< mpV(1)<<", "<< mpV(2)<<") "<< endl
1579 << "n2 = m.AddNode( "<< mp[0](0)<<", "<< mp[0](1)<<", "<< mp[0](2)<<") "<< endl
1580 << "n3 = m.AddNode( "<< mp[1](0)<<", "<< mp[1](1)<<", "<< mp[1](2)<<" )" << endl
1581 << "m.AddFace([n1,n2,n3])" << endl;
1583 } // loop on internal vertices of a solid
1585 } // loop on solids with internal vertices
1588 //================================================================================
1590 * \brief Fill netgen mesh with segments of a FACE
1591 * \param ngMesh - netgen mesh
1592 * \param geom - container of OCCT geometry to mesh
1593 * \param wires - data of nodes on FACE boundary
1594 * \param helper - mesher helper holding the FACE
1595 * \param nodeVec - vector of nodes in which node index == netgen ID
1596 * \retval SMESH_ComputeErrorPtr - error description
1598 //================================================================================
1600 SMESH_ComputeErrorPtr
1601 NETGENPlugin_Mesher::AddSegmentsToMesh(netgen::Mesh& ngMesh,
1602 netgen::OCCGeometry& geom,
1603 const TSideVector& wires,
1604 SMESH_MesherHelper& helper,
1605 vector< const SMDS_MeshNode* > & nodeVec,
1606 const bool overrideMinH)
1608 // ----------------------------
1609 // Check wires and count nodes
1610 // ----------------------------
1612 for ( int iW = 0; iW < wires.size(); ++iW )
1614 StdMeshers_FaceSidePtr wire = wires[ iW ];
1615 if ( wire->MissVertexNode() )
1617 // Commented for issue 0020960. It worked for the case, let's wait for case where it doesn't.
1618 // It seems that there is no reason for this limitation
1620 // (new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH, "Missing nodes on vertices"));
1622 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1623 if ( uvPtVec.size() != wire->NbPoints() )
1624 return SMESH_ComputeError::New(COMPERR_BAD_INPUT_MESH,
1625 SMESH_Comment("Unexpected nb of points on wire ") << iW
1626 << ": " << uvPtVec.size()<<" != "<<wire->NbPoints());
1627 nbNodes += wire->NbPoints();
1629 nodeVec.reserve( nodeVec.size() + nbNodes + 1 );
1630 if ( nodeVec.empty() )
1631 nodeVec.push_back( 0 );
1633 // -----------------
1635 // -----------------
1637 const bool wasNgMeshEmpty = ( ngMesh.GetNP() < 1 ); /* true => this method is called by
1638 NETGENPlugin_NETGEN_2D_ONLY */
1640 // map for nodes on vertices since they can be shared between wires
1641 // ( issue 0020676, face_int_box.brep) and nodes built by NETGEN
1642 map<const SMDS_MeshNode*, int > node2ngID;
1643 if ( !wasNgMeshEmpty ) // fill node2ngID with nodes built by NETGEN
1645 set< int > subIDs; // ids of sub-shapes of the FACE
1646 for ( int iW = 0; iW < wires.size(); ++iW )
1648 StdMeshers_FaceSidePtr wire = wires[ iW ];
1649 for ( int iE = 0, nbE = wire->NbEdges(); iE < nbE; ++iE )
1651 subIDs.insert( wire->EdgeID( iE ));
1652 subIDs.insert( helper.GetMeshDS()->ShapeToIndex( wire->FirstVertex( iE )));
1655 for ( size_t ngID = 1; ngID < nodeVec.size(); ++ngID )
1656 if ( subIDs.count( nodeVec[ngID]->getshapeId() ))
1657 node2ngID.insert( make_pair( nodeVec[ngID], ngID ));
1660 const int solidID = 0, faceID = geom.fmap.FindIndex( helper.GetSubShape() );
1661 if ( ngMesh.GetNFD() < 1 )
1662 ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(faceID, solidID, solidID, 0));
1664 for ( int iW = 0; iW < wires.size(); ++iW )
1666 StdMeshers_FaceSidePtr wire = wires[ iW ];
1667 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1668 const int nbSegments = wire->NbPoints() - 1;
1670 // assure the 1st node to be in node2ngID, which is needed to correctly
1671 // "close chain of segments" (see below) in case if the 1st node is not
1672 // onVertex because it is on a Viscous layer
1673 node2ngID.insert( make_pair( uvPtVec[ 0 ].node, ngMesh.GetNP() + 1 ));
1675 // compute length of every segment
1676 vector<double> segLen( nbSegments );
1677 for ( int i = 0; i < nbSegments; ++i )
1678 segLen[i] = SMESH_TNodeXYZ( uvPtVec[ i ].node ).Distance( uvPtVec[ i+1 ].node );
1680 int edgeID = 1, posID = -2;
1681 bool isInternalWire = false;
1682 double vertexNormPar = 0;
1683 const int prevNbNGSeg = ngMesh.GetNSeg();
1684 for ( int i = 0; i < nbSegments; ++i ) // loop on segments
1686 // Add the first point of a segment
1688 const SMDS_MeshNode * n = uvPtVec[ i ].node;
1689 const int posShapeID = n->getshapeId();
1690 bool onVertex = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX );
1691 bool onEdge = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_EDGE );
1693 // skip nodes on degenerated edges
1694 if ( helper.IsDegenShape( posShapeID ) &&
1695 helper.IsDegenShape( uvPtVec[ i+1 ].node->getshapeId() ))
1698 int ngID1 = ngMesh.GetNP() + 1, ngID2 = ngID1+1;
1699 if ( onVertex || ( !wasNgMeshEmpty && onEdge ))
1700 ngID1 = node2ngID.insert( make_pair( n, ngID1 )).first->second;
1701 if ( ngID1 > ngMesh.GetNP() )
1703 netgen::MeshPoint mp( netgen::Point<3> (n->X(), n->Y(), n->Z()) );
1704 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1705 nodeVec.push_back( n );
1707 else // n is in ngMesh already, and ngID2 in prev segment is wrong
1709 ngID2 = ngMesh.GetNP() + 1;
1710 if ( i > 0 ) // prev segment belongs to same wire
1712 netgen::Segment& prevSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
1719 netgen::Segment seg;
1721 seg[0] = ngID1; // ng node id
1722 seg[1] = ngID2; // ng node id
1723 seg.edgenr = ngMesh.GetNSeg() + 1; // ng segment id
1724 seg.si = faceID; // = geom.fmap.FindIndex (face);
1726 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1728 const UVPtStruct& pnt = uvPtVec[ i + iEnd ];
1730 seg.epgeominfo[ iEnd ].dist = pnt.param; // param on curve
1731 seg.epgeominfo[ iEnd ].u = pnt.u;
1732 seg.epgeominfo[ iEnd ].v = pnt.v;
1734 // find out edge id and node parameter on edge
1735 onVertex = ( pnt.normParam + 1e-10 > vertexNormPar );
1736 if ( onVertex || posShapeID != posID )
1739 double normParam = pnt.normParam;
1741 normParam = 0.5 * ( uvPtVec[ i ].normParam + uvPtVec[ i+1 ].normParam );
1742 int edgeIndexInWire = wire->EdgeIndex( normParam );
1743 vertexNormPar = wire->LastParameter( edgeIndexInWire );
1744 const TopoDS_Edge& edge = wire->Edge( edgeIndexInWire );
1745 edgeID = geom.emap.FindIndex( edge );
1747 isInternalWire = ( edge.Orientation() == TopAbs_INTERNAL );
1748 // if ( onVertex ) // param on curve is different on each of two edges
1749 // seg.epgeominfo[ iEnd ].dist = helper.GetNodeU( edge, pnt.node );
1751 seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
1754 ngMesh.AddSegment (seg);
1756 // restrict size of elements near the segment
1757 SMESH_TNodeXYZ np1( n ), np2( uvPtVec[ i+1 ].node );
1758 // get an average size of adjacent segments to avoid sharp change of
1759 // element size (regression on issue 0020452, note 0010898)
1760 int iPrev = SMESH_MesherHelper::WrapIndex( i-1, nbSegments );
1761 int iNext = SMESH_MesherHelper::WrapIndex( i+1, nbSegments );
1762 double sumH = segLen[ iPrev ] + segLen[ i ] + segLen[ iNext ];
1763 int nbSeg = ( int( segLen[ iPrev ] > sumH / 100.) +
1764 int( segLen[ i ] > sumH / 100.) +
1765 int( segLen[ iNext ] > sumH / 100.));
1767 RestrictLocalSize( ngMesh, 0.5*(np1+np2), sumH / nbSeg, overrideMinH );
1769 if ( isInternalWire )
1771 swap (seg[0], seg[1]);
1772 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1773 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1774 ngMesh.AddSegment (seg);
1776 } // loop on segments on a wire
1778 // close chain of segments
1779 if ( nbSegments > 0 )
1781 netgen::Segment& lastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() - int( isInternalWire));
1782 const SMDS_MeshNode * lastNode = uvPtVec.back().node;
1783 lastSeg[1] = node2ngID.insert( make_pair( lastNode, lastSeg[1] )).first->second;
1784 if ( lastSeg[1] > ngMesh.GetNP() )
1786 netgen::MeshPoint mp( netgen::Point<3> (lastNode->X(), lastNode->Y(), lastNode->Z()) );
1787 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1788 nodeVec.push_back( lastNode );
1790 if ( isInternalWire )
1792 netgen::Segment& realLastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
1793 realLastSeg[0] = lastSeg[1];
1797 #ifdef DUMP_SEGMENTS
1798 cout << "BEGIN WIRE " << iW << endl;
1799 for ( int i = prevNbNGSeg+1; i <= ngMesh.GetNSeg(); ++i )
1801 netgen::Segment& seg = ngMesh.LineSegment( i );
1803 netgen::Segment& prevSeg = ngMesh.LineSegment( i-1 );
1804 if ( seg[0] == prevSeg[1] && seg[1] == prevSeg[0] )
1806 cout << "Segment: " << seg.edgenr << endl << "\tis REVRESE of the previous one" << endl;
1810 cout << "Segment: " << seg.edgenr << endl
1811 << "\tp1: " << seg[0] << endl
1812 << "\tp2: " << seg[1] << endl
1813 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
1814 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
1815 << "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
1816 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
1817 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl
1818 << "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
1820 cout << "--END WIRE " << iW << endl;
1823 } // loop on WIREs of a FACE
1825 // add a segment instead of an internal vertex
1826 if ( wasNgMeshEmpty )
1828 NETGENPlugin_Internals intShapes( *helper.GetMesh(), helper.GetSubShape(), /*is3D=*/false );
1829 AddIntVerticesInFaces( geom, ngMesh, nodeVec, intShapes );
1831 ngMesh.CalcSurfacesOfNode();
1836 //================================================================================
1838 * \brief Fill SMESH mesh according to contents of netgen mesh
1839 * \param occgeo - container of OCCT geometry to mesh
1840 * \param ngMesh - netgen mesh
1841 * \param initState - bn of entities in netgen mesh before computing
1842 * \param sMesh - SMESH mesh to fill in
1843 * \param nodeVec - vector of nodes in which node index == netgen ID
1844 * \param comment - returns problem description
1845 * \param quadHelper - holder of medium nodes of sub-meshes
1846 * \retval int - error
1848 //================================================================================
1850 int NETGENPlugin_Mesher::FillSMesh(const netgen::OCCGeometry& occgeo,
1851 netgen::Mesh& ngMesh,
1852 const NETGENPlugin_ngMeshInfo& initState,
1854 std::vector<const SMDS_MeshNode*>& nodeVec,
1855 SMESH_Comment& comment,
1856 SMESH_MesherHelper* quadHelper)
1858 int nbNod = ngMesh.GetNP();
1859 int nbSeg = ngMesh.GetNSeg();
1860 int nbFac = ngMesh.GetNSE();
1861 int nbVol = ngMesh.GetNE();
1863 SMESHDS_Mesh* meshDS = sMesh.GetMeshDS();
1865 // quadHelper is used for either
1866 // 1) making quadratic elements when a lower dimention mesh is loaded
1867 // to SMESH before convertion to quadratic by NETGEN
1868 // 2) sewing of quadratic elements with quadratic elements of sub-meshes
1869 if ( quadHelper && !quadHelper->GetIsQuadratic() && quadHelper->GetTLinkNodeMap().empty() )
1872 // -------------------------------------
1873 // Create and insert nodes into nodeVec
1874 // -------------------------------------
1876 nodeVec.resize( nbNod + 1 );
1877 int i, nbInitNod = initState._nbNodes;
1878 for (i = nbInitNod+1; i <= nbNod; ++i )
1880 const netgen::MeshPoint& ngPoint = ngMesh.Point(i);
1881 SMDS_MeshNode* node = NULL;
1882 TopoDS_Vertex aVert;
1883 // First, netgen creates nodes on vertices in occgeo.vmap,
1884 // so node index corresponds to vertex index
1885 // but (issue 0020776) netgen does not create nodes with equal coordinates
1886 if ( i-nbInitNod <= occgeo.vmap.Extent() )
1888 gp_Pnt p ( NGPOINT_COORDS(ngPoint) );
1889 for (int iV = i-nbInitNod; aVert.IsNull() && iV <= occgeo.vmap.Extent(); ++iV)
1891 aVert = TopoDS::Vertex( occgeo.vmap( iV ));
1892 gp_Pnt pV = BRep_Tool::Pnt( aVert );
1893 if ( p.SquareDistance( pV ) > 1e-20 )
1896 node = const_cast<SMDS_MeshNode*>( SMESH_Algo::VertexNode( aVert, meshDS ));
1899 if (!node) // node not found on vertex
1901 node = meshDS->AddNode( NGPOINT_COORDS( ngPoint ));
1902 if (!aVert.IsNull())
1903 meshDS->SetNodeOnVertex(node, aVert);
1908 // -------------------------------------------
1909 // Create mesh segments along geometric edges
1910 // -------------------------------------------
1912 int nbInitSeg = initState._nbSegments;
1913 for (i = nbInitSeg+1; i <= nbSeg; ++i )
1915 const netgen::Segment& seg = ngMesh.LineSegment(i);
1917 int pinds[3] = { seg.pnums[0], seg.pnums[1], seg.pnums[2] };
1920 for (int j=0; j < 3; ++j)
1922 int pind = pinds[j];
1923 if (pind <= 0 || !nodeVec_ACCESS(pind))
1931 int aGeomEdgeInd = seg.epgeominfo[j].edgenr;
1932 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
1933 aEdge = TopoDS::Edge(occgeo.emap(aGeomEdgeInd));
1935 param = seg.epgeominfo[j].dist;
1938 else // middle point
1940 param = param2 * 0.5;
1942 if (!aEdge.IsNull() && nodeVec_ACCESS(pind)->getshapeId() < 1)
1944 meshDS->SetNodeOnEdge(nodeVec_ACCESS(pind), aEdge, param);
1949 SMDS_MeshEdge* edge = 0;
1950 if (nbp == 2) // second order ?
1952 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1])))
1954 if ( quadHelper ) // final mesh must be quadratic
1955 edge = quadHelper->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
1957 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
1961 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
1962 nodeVec_ACCESS(pinds[2])))
1964 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
1965 nodeVec_ACCESS(pinds[2]));
1969 if ( comment.empty() ) comment << "Cannot create a mesh edge";
1970 MESSAGE("Cannot create a mesh edge");
1971 nbSeg = nbFac = nbVol = 0;
1974 if ( !aEdge.IsNull() && edge->getshapeId() < 1 )
1975 meshDS->SetMeshElementOnShape(edge, aEdge);
1977 else if ( comment.empty() )
1979 comment << "Invalid netgen segment #" << i;
1983 // ----------------------------------------
1984 // Create mesh faces along geometric faces
1985 // ----------------------------------------
1987 int nbInitFac = initState._nbFaces;
1988 int quadFaceID = ngMesh.GetNFD() + 1;
1989 if ( nbInitFac < nbFac )
1990 // add a faces descriptor to exclude qudrangle elements generated by NETGEN
1991 // from computation of 3D mesh
1992 ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(quadFaceID, /*solid1=*/0, /*solid2=*/0, 0));
1994 vector<const SMDS_MeshNode*> nodes;
1995 for (i = nbInitFac+1; i <= nbFac; ++i )
1997 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1998 int aGeomFaceInd = elem.GetIndex();
2000 if (aGeomFaceInd > 0 && aGeomFaceInd <= occgeo.fmap.Extent())
2001 aFace = TopoDS::Face(occgeo.fmap(aGeomFaceInd));
2003 for (int j=1; j <= elem.GetNP(); ++j)
2005 int pind = elem.PNum(j);
2006 if ( pind < 1 || pind >= nodeVec.size() )
2008 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind))
2010 nodes.push_back( node );
2011 if (!aFace.IsNull() && node->getshapeId() < 1)
2013 const netgen::PointGeomInfo& pgi = elem.GeomInfoPi(j);
2014 meshDS->SetNodeOnFace(node, aFace, pgi.u, pgi.v);
2018 if ( nodes.size() != elem.GetNP() )
2020 if ( comment.empty() )
2021 comment << "Invalid netgen 2d element #" << i;
2022 continue; // bad node ids
2024 SMDS_MeshFace* face = NULL;
2025 switch (elem.GetType())
2028 if ( quadHelper ) // final mesh must be quadratic
2029 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2]);
2031 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2]);
2034 if ( quadHelper ) // final mesh must be quadratic
2035 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2037 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2038 // exclude qudrangle elements from computation of 3D mesh
2039 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2042 nodes[5] = mediumNode( nodes[0],nodes[1],nodes[5], quadHelper );
2043 nodes[3] = mediumNode( nodes[1],nodes[2],nodes[3], quadHelper );
2044 nodes[4] = mediumNode( nodes[2],nodes[0],nodes[4], quadHelper );
2045 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[5],nodes[3],nodes[4]);
2048 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2049 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2050 nodes[5] = mediumNode( nodes[2],nodes[3],nodes[5], quadHelper );
2051 nodes[6] = mediumNode( nodes[3],nodes[0],nodes[6], quadHelper );
2052 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3],
2053 nodes[4],nodes[7],nodes[5],nodes[6]);
2054 // exclude qudrangle elements from computation of 3D mesh
2055 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2058 MESSAGE("NETGEN created a face of unexpected type, ignoring");
2063 if ( comment.empty() ) comment << "Cannot create a mesh face";
2064 MESSAGE("Cannot create a mesh face");
2065 nbSeg = nbFac = nbVol = 0;
2068 if (!aFace.IsNull())
2069 meshDS->SetMeshElementOnShape(face, aFace);
2072 // ------------------
2073 // Create tetrahedra
2074 // ------------------
2076 for (i = 1; i <= nbVol; ++i)
2078 const netgen::Element& elem = ngMesh.VolumeElement(i);
2079 int aSolidInd = elem.GetIndex();
2080 TopoDS_Solid aSolid;
2081 if (aSolidInd > 0 && aSolidInd <= occgeo.somap.Extent())
2082 aSolid = TopoDS::Solid(occgeo.somap(aSolidInd));
2084 for (int j=1; j <= elem.GetNP(); ++j)
2086 int pind = elem.PNum(j);
2087 if ( pind < 1 || pind >= nodeVec.size() )
2089 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind) )
2091 nodes.push_back(node);
2092 if ( !aSolid.IsNull() && node->getshapeId() < 1 )
2093 meshDS->SetNodeInVolume(node, aSolid);
2096 if ( nodes.size() != elem.GetNP() )
2098 if ( comment.empty() )
2099 comment << "Invalid netgen 3d element #" << i;
2102 SMDS_MeshVolume* vol = NULL;
2103 switch (elem.GetType())
2106 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3]);
2109 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2110 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2111 nodes[5] = mediumNode( nodes[2],nodes[0],nodes[5], quadHelper );
2112 nodes[6] = mediumNode( nodes[0],nodes[3],nodes[6], quadHelper );
2113 nodes[8] = mediumNode( nodes[1],nodes[3],nodes[8], quadHelper );
2114 nodes[9] = mediumNode( nodes[2],nodes[3],nodes[9], quadHelper );
2115 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3],
2116 nodes[4],nodes[7],nodes[5],nodes[6],nodes[8],nodes[9]);
2119 MESSAGE("NETGEN created a volume of unexpected type, ignoring");
2124 if ( comment.empty() ) comment << "Cannot create a mesh volume";
2125 MESSAGE("Cannot create a mesh volume");
2126 nbSeg = nbFac = nbVol = 0;
2129 if (!aSolid.IsNull())
2130 meshDS->SetMeshElementOnShape(vol, aSolid);
2132 return comment.empty() ? 0 : 1;
2137 //================================================================================
2139 * \brief Restrict size of elements on the given edge
2141 //================================================================================
2143 void setLocalSize(const TopoDS_Edge& edge,
2147 const int nb = 1000;
2148 Standard_Real u1, u2;
2149 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, u1, u2);
2150 if ( curve.IsNull() )
2152 TopoDS_Iterator vIt( edge );
2153 if ( !vIt.More() ) return;
2154 gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex( vIt.Value() ));
2155 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
2159 Standard_Real delta = (u2-u1)/nb;
2160 for(int i=0; i<nb; i++)
2162 Standard_Real u = u1 + delta*i;
2163 gp_Pnt p = curve->Value(u);
2164 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
2165 netgen::Point3d pi(p.X(), p.Y(), p.Z());
2166 double resultSize = mesh.GetH(pi);
2167 if ( resultSize - size > 0.1*size )
2168 // netgen does restriction iff oldH/newH > 1.2 (localh.cpp:136)
2169 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), resultSize/1.201 );
2174 //================================================================================
2176 * \brief Convert error into text
2178 //================================================================================
2180 std::string text(int err)
2185 SMESH_Comment("Error in netgen::OCCGenerateMesh() at ") << netgen::multithread.task;
2188 //================================================================================
2190 * \brief Convert exception into text
2192 //================================================================================
2194 std::string text(Standard_Failure& ex)
2196 SMESH_Comment str("Exception in netgen::OCCGenerateMesh()");
2197 str << " at " << netgen::multithread.task
2198 << ": " << ex.DynamicType()->Name();
2199 if ( ex.GetMessageString() && strlen( ex.GetMessageString() ))
2200 str << ": " << ex.GetMessageString();
2203 //================================================================================
2205 * \brief Convert exception into text
2207 //================================================================================
2209 std::string text(netgen::NgException& ex)
2211 SMESH_Comment str("NgException");
2212 if ( strlen( netgen::multithread.task ) > 0 )
2213 str << " at " << netgen::multithread.task;
2214 str << ": " << ex.What();
2218 //================================================================================
2220 * \brief Looks for triangles lying on a SOLID
2222 //================================================================================
2224 bool hasBadElemOnSolid( const list<const SMDS_MeshElement*>& elems,
2225 SMESH_subMesh* solidSM )
2227 TopTools_IndexedMapOfShape solidSubs;
2228 TopExp::MapShapes( solidSM->GetSubShape(), solidSubs );
2229 SMESHDS_Mesh* mesh = solidSM->GetFather()->GetMeshDS();
2231 list<const SMDS_MeshElement*>::const_iterator e = elems.begin();
2232 for ( ; e != elems.end(); ++e )
2234 const SMDS_MeshElement* elem = *e;
2235 // if ( elem->GetType() != SMDSAbs_Face ) -- 23047
2237 int nbNodesOnSolid = 0, nbNodes = elem->NbNodes();
2238 SMDS_NodeIteratorPtr nIt = elem->nodeIterator();
2239 while ( nIt->more() )
2241 const SMDS_MeshNode* n = nIt->next();
2242 const TopoDS_Shape& s = mesh->IndexToShape( n->getshapeId() );
2243 nbNodesOnSolid += ( !s.IsNull() && solidSubs.Contains( s ));
2244 if ( nbNodesOnSolid > 2 ||
2245 nbNodesOnSolid == nbNodes)
2252 const double edgeMeshingTime = 0.001;
2253 const double faceMeshingTime = 0.019;
2254 const double edgeFaceMeshingTime = edgeMeshingTime + faceMeshingTime;
2255 const double faceOptimizTime = 0.06;
2256 const double voluMeshingTime = 0.15;
2257 const double volOptimizeTime = 0.77;
2260 //=============================================================================
2262 * Here we are going to use the NETGEN mesher
2264 //=============================================================================
2266 bool NETGENPlugin_Mesher::Compute()
2268 NETGENPlugin_NetgenLibWrapper ngLib;
2270 netgen::MeshingParameters& mparams = netgen::mparam;
2271 MESSAGE("Compute with:\n"
2272 " max size = " << mparams.maxh << "\n"
2273 " segments per edge = " << mparams.segmentsperedge);
2275 " growth rate = " << mparams.grading << "\n"
2276 " elements per radius = " << mparams.curvaturesafety << "\n"
2277 " second order = " << mparams.secondorder << "\n"
2278 " quad allowed = " << mparams.quad << "\n"
2279 " surface curvature = " << mparams.uselocalh << "\n"
2280 " fuse edges = " << netgen::merge_solids);
2282 SMESH_ComputeErrorPtr error = SMESH_ComputeError::New();
2283 SMESH_MesherHelper quadHelper( *_mesh );
2284 quadHelper.SetIsQuadratic( mparams.secondorder );
2286 static string debugFile = "/tmp/ngMesh.py"; /* to call toPython( _ngMesh, debugFile )
2287 while debugging netgen */
2288 // -------------------------
2289 // Prepare OCC geometry
2290 // -------------------------
2292 netgen::OCCGeometry occgeo;
2293 list< SMESH_subMesh* > meshedSM[3]; // for 0-2 dimensions
2294 NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
2295 PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
2298 _totalTime = edgeFaceMeshingTime;
2300 _totalTime += faceOptimizTime;
2302 _totalTime += voluMeshingTime + ( _optimize ? volOptimizeTime : 0 );
2303 double doneTime = 0;
2306 _curShapeIndex = -1;
2308 // -------------------------
2309 // Generate the mesh
2310 // -------------------------
2313 NETGENPlugin_ngMeshInfo initState; // it remembers size of ng mesh equal to size of Smesh
2315 SMESH_Comment comment;
2318 // vector of nodes in which node index == netgen ID
2319 vector< const SMDS_MeshNode* > nodeVec;
2327 // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
2328 mparams.uselocalh = false;
2329 mparams.grading = 0.8; // not limitited size growth
2331 if ( _simpleHyp->GetNumberOfSegments() )
2333 mparams.maxh = occgeo.boundingbox.Diam();
2336 mparams.maxh = _simpleHyp->GetLocalLength();
2339 if ( mparams.maxh == 0.0 )
2340 mparams.maxh = occgeo.boundingbox.Diam();
2341 if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
2342 mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
2344 // Local size on faces
2345 occgeo.face_maxh = mparams.maxh;
2347 // Let netgen create _ngMesh and calculate element size on not meshed shapes
2351 int startWith = netgen::MESHCONST_ANALYSE;
2352 int endWith = netgen::MESHCONST_ANALYSE;
2357 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2359 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2361 if(netgen::multithread.terminate)
2364 comment << text(err);
2366 catch (Standard_Failure& ex)
2368 comment << text(ex);
2370 err = 0; //- MESHCONST_ANALYSE isn't so important step
2373 ngLib.setMesh(( Ng_Mesh*) _ngMesh );
2375 _ngMesh->ClearFaceDescriptors(); // we make descriptors our-self
2379 // Pass 1D simple parameters to NETGEN
2380 // --------------------------------
2381 int nbSeg = _simpleHyp->GetNumberOfSegments();
2382 double segSize = _simpleHyp->GetLocalLength();
2383 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
2385 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
2387 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
2388 setLocalSize( e, segSize, *_ngMesh );
2391 else // if ( ! _simpleHyp )
2393 // Local size on vertices and edges
2394 // --------------------------------
2395 for(std::map<int,double>::const_iterator it=EdgeId2LocalSize.begin(); it!=EdgeId2LocalSize.end(); it++)
2397 int key = (*it).first;
2398 double hi = (*it).second;
2399 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2400 const TopoDS_Edge& e = TopoDS::Edge(shape);
2401 setLocalSize( e, hi, *_ngMesh );
2403 for(std::map<int,double>::const_iterator it=VertexId2LocalSize.begin(); it!=VertexId2LocalSize.end(); it++)
2405 int key = (*it).first;
2406 double hi = (*it).second;
2407 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2408 const TopoDS_Vertex& v = TopoDS::Vertex(shape);
2409 gp_Pnt p = BRep_Tool::Pnt(v);
2410 NETGENPlugin_Mesher::RestrictLocalSize( *_ngMesh, p.XYZ(), hi );
2412 for(map<int,double>::const_iterator it=FaceId2LocalSize.begin();
2413 it!=FaceId2LocalSize.end(); it++)
2415 int key = (*it).first;
2416 double val = (*it).second;
2417 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2418 int faceNgID = occgeo.fmap.FindIndex(shape);
2419 occgeo.SetFaceMaxH(faceNgID, val);
2420 for ( TopExp_Explorer edgeExp( shape, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
2421 setLocalSize( TopoDS::Edge( edgeExp.Current() ), val, *_ngMesh );
2425 // Precompute internal edges (issue 0020676) in order to
2426 // add mesh on them correctly (twice) to netgen mesh
2427 if ( !err && internals.hasInternalEdges() )
2429 // load internal shapes into OCCGeometry
2430 netgen::OCCGeometry intOccgeo;
2431 internals.getInternalEdges( intOccgeo.fmap, intOccgeo.emap, intOccgeo.vmap, meshedSM );
2432 intOccgeo.boundingbox = occgeo.boundingbox;
2433 intOccgeo.shape = occgeo.shape;
2434 intOccgeo.face_maxh.SetSize(intOccgeo.fmap.Extent());
2435 intOccgeo.face_maxh = netgen::mparam.maxh;
2436 netgen::Mesh *tmpNgMesh = NULL;
2440 // compute local H on internal shapes in the main mesh
2441 //OCCSetLocalMeshSize(intOccgeo, *_ngMesh); it deletes _ngMesh->localH
2443 // let netgen create a temporary mesh
2445 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2447 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2449 if(netgen::multithread.terminate)
2452 // copy LocalH from the main to temporary mesh
2453 initState.transferLocalH( _ngMesh, tmpNgMesh );
2455 // compute mesh on internal edges
2456 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2458 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2460 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2462 comment << text(err);
2464 catch (Standard_Failure& ex)
2466 comment << text(ex);
2469 initState.restoreLocalH( tmpNgMesh );
2471 // fill SMESH by netgen mesh
2472 vector< const SMDS_MeshNode* > tmpNodeVec;
2473 FillSMesh( intOccgeo, *tmpNgMesh, initState, *_mesh, tmpNodeVec, comment );
2474 err = ( err || !comment.empty() );
2476 nglib::Ng_DeleteMesh((nglib::Ng_Mesh*)tmpNgMesh);
2479 // Fill _ngMesh with nodes and segments of computed submeshes
2482 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_0D ]) &&
2483 FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_1D ], &quadHelper));
2485 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2490 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2495 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2497 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2499 if(netgen::multithread.terminate)
2502 comment << text(err);
2504 catch (Standard_Failure& ex)
2506 comment << text(ex);
2511 _ticTime = ( doneTime += edgeMeshingTime ) / _totalTime / _progressTic;
2513 mparams.uselocalh = true; // restore as it is used at surface optimization
2515 // ---------------------
2516 // compute surface mesh
2517 // ---------------------
2520 // Pass 2D simple parameters to NETGEN
2522 if ( double area = _simpleHyp->GetMaxElementArea() ) {
2524 mparams.maxh = sqrt(2. * area/sqrt(3.0));
2525 mparams.grading = 0.4; // moderate size growth
2528 // length from edges
2529 if ( _ngMesh->GetNSeg() ) {
2530 double edgeLength = 0;
2531 TopTools_MapOfShape visitedEdges;
2532 for ( TopExp_Explorer exp( _shape, TopAbs_EDGE ); exp.More(); exp.Next() )
2533 if( visitedEdges.Add(exp.Current()) )
2534 edgeLength += SMESH_Algo::EdgeLength( TopoDS::Edge( exp.Current() ));
2535 // we have to multiply length by 2 since for each TopoDS_Edge there
2536 // are double set of NETGEN edges, in other words, we have to
2537 // divide _ngMesh->GetNSeg() by 2.
2538 mparams.maxh = 2*edgeLength / _ngMesh->GetNSeg();
2541 mparams.maxh = 1000;
2543 mparams.grading = 0.2; // slow size growth
2545 mparams.quad = _simpleHyp->GetAllowQuadrangles();
2546 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2547 _ngMesh->SetGlobalH (mparams.maxh);
2548 netgen::Box<3> bb = occgeo.GetBoundingBox();
2549 bb.Increase (bb.Diam()/20);
2550 _ngMesh->SetLocalH (bb.PMin(), bb.PMax(), mparams.grading);
2553 // Care of vertices internal in faces (issue 0020676)
2554 if ( internals.hasInternalVertexInFace() )
2556 // store computed segments in SMESH in order not to create SMESH
2557 // edges for ng segments added by AddIntVerticesInFaces()
2558 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2559 // add segments to faces with internal vertices
2560 AddIntVerticesInFaces( occgeo, *_ngMesh, nodeVec, internals );
2561 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2564 // Build viscous layers
2565 if ( _isViscousLayers2D )
2567 if ( !internals.hasInternalVertexInFace() ) {
2568 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2569 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2571 SMESH_ProxyMesh::Ptr viscousMesh;
2572 SMESH_MesherHelper helper( *_mesh );
2573 for ( int faceID = 1; faceID <= occgeo.fmap.Extent(); ++faceID )
2575 const TopoDS_Face& F = TopoDS::Face( occgeo.fmap( faceID ));
2576 viscousMesh = StdMeshers_ViscousLayers2D::Compute( *_mesh, F );
2579 // exclude from computation ng segments built on EDGEs of F
2580 for (int i = 1; i <= _ngMesh->GetNSeg(); i++)
2582 netgen::Segment & seg = _ngMesh->LineSegment(i);
2583 if (seg.si == faceID)
2586 // add new segments to _ngMesh instead of excluded ones
2587 helper.SetSubShape( F );
2589 StdMeshers_FaceSide::GetFaceWires( F, *_mesh, /*skipMediumNodes=*/true,
2590 error, viscousMesh );
2591 error = AddSegmentsToMesh( *_ngMesh, occgeo, wires, helper, nodeVec );
2593 if ( !error ) error = SMESH_ComputeError::New();
2595 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2598 // Let netgen compute 2D mesh
2599 startWith = netgen::MESHCONST_MESHSURFACE;
2600 endWith = _optimize ? netgen::MESHCONST_OPTSURFACE : netgen::MESHCONST_MESHSURFACE;
2605 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2607 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2609 if(netgen::multithread.terminate)
2612 comment << text (err);
2614 catch (Standard_Failure& ex)
2616 comment << text(ex);
2617 //err = 1; -- try to make volumes anyway
2619 catch (netgen::NgException exc)
2621 comment << text(exc);
2622 //err = 1; -- try to make volumes anyway
2627 doneTime += faceMeshingTime + ( _optimize ? faceOptimizTime : 0 );
2628 _ticTime = doneTime / _totalTime / _progressTic;
2630 // ---------------------
2631 // generate volume mesh
2632 // ---------------------
2633 // Fill _ngMesh with nodes and faces of computed 2D submeshes
2634 if ( !err && _isVolume && ( !meshedSM[ MeshDim_2D ].empty() || mparams.quad ))
2636 // load SMESH with computed segments and faces
2637 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2639 // compute pyramids on quadrangles
2640 SMESH_ProxyMesh::Ptr proxyMesh;
2641 if ( _mesh->NbQuadrangles() > 0 )
2642 for ( int iS = 1; iS <= occgeo.somap.Extent(); ++iS )
2644 StdMeshers_QuadToTriaAdaptor* Adaptor = new StdMeshers_QuadToTriaAdaptor;
2645 proxyMesh.reset( Adaptor );
2647 int nbPyrams = _mesh->NbPyramids();
2648 Adaptor->Compute( *_mesh, occgeo.somap(iS) );
2649 if ( nbPyrams != _mesh->NbPyramids() )
2651 list< SMESH_subMesh* > quadFaceSM;
2652 for (TopExp_Explorer face(occgeo.somap(iS), TopAbs_FACE); face.More(); face.Next())
2653 if ( Adaptor->GetProxySubMesh( face.Current() ))
2655 quadFaceSM.push_back( _mesh->GetSubMesh( face.Current() ));
2656 meshedSM[ MeshDim_2D ].remove( quadFaceSM.back() );
2658 FillNgMesh(occgeo, *_ngMesh, nodeVec, quadFaceSM, &quadHelper, proxyMesh);
2661 // fill _ngMesh with faces of sub-meshes
2662 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_2D ], &quadHelper));
2663 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2664 //toPython( _ngMesh, "/tmp/ngPython.py");
2666 if (!err && _isVolume)
2668 // Pass 3D simple parameters to NETGEN
2669 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
2670 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
2672 if ( double vol = simple3d->GetMaxElementVolume() ) {
2674 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
2675 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2678 // length from faces
2679 mparams.maxh = _ngMesh->AverageH();
2681 _ngMesh->SetGlobalH (mparams.maxh);
2682 mparams.grading = 0.4;
2684 _ngMesh->CalcLocalH(mparams.grading);
2686 _ngMesh->CalcLocalH();
2689 // Care of vertices internal in solids and internal faces (issue 0020676)
2690 if ( internals.hasInternalVertexInSolid() || internals.hasInternalFaces() )
2692 // store computed faces in SMESH in order not to create SMESH
2693 // faces for ng faces added here
2694 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2695 // add ng faces to solids with internal vertices
2696 AddIntVerticesInSolids( occgeo, *_ngMesh, nodeVec, internals );
2697 // duplicate mesh faces on internal faces
2698 FixIntFaces( occgeo, *_ngMesh, internals );
2699 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2701 // Let netgen compute 3D mesh
2702 startWith = endWith = netgen::MESHCONST_MESHVOLUME;
2707 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2709 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2711 if(netgen::multithread.terminate)
2714 if ( comment.empty() ) // do not overwrite a previos error
2715 comment << text(err);
2717 catch (Standard_Failure& ex)
2719 if ( comment.empty() ) // do not overwrite a previos error
2720 comment << text(ex);
2723 catch (netgen::NgException exc)
2725 if ( comment.empty() ) // do not overwrite a previos error
2726 comment << text(exc);
2729 _ticTime = ( doneTime += voluMeshingTime ) / _totalTime / _progressTic;
2731 // Let netgen optimize 3D mesh
2732 if ( !err && _optimize )
2734 startWith = endWith = netgen::MESHCONST_OPTVOLUME;
2739 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2741 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2743 if(netgen::multithread.terminate)
2746 if ( comment.empty() ) // do not overwrite a previos error
2747 comment << text(err);
2749 catch (Standard_Failure& ex)
2751 if ( comment.empty() ) // do not overwrite a previos error
2752 comment << text(ex);
2754 catch (netgen::NgException exc)
2756 if ( comment.empty() ) // do not overwrite a previos error
2757 comment << text(exc);
2761 if (!err && mparams.secondorder > 0)
2766 if ( !meshedSM[ MeshDim_1D ].empty() )
2768 // remove segments not attached to geometry (IPAL0052479)
2769 for (int i = 1; i <= _ngMesh->GetNSeg(); ++i)
2771 const netgen::Segment & seg = _ngMesh->LineSegment (i);
2772 if ( seg.epgeominfo[ 0 ].edgenr == 0 )
2773 _ngMesh->DeleteSegment( i );
2775 _ngMesh->Compress();
2777 // convert to quadratic
2778 netgen::OCCRefinementSurfaces ref (occgeo);
2779 ref.MakeSecondOrder (*_ngMesh);
2781 // care of elements already loaded to SMESH
2782 // if ( initState._nbSegments > 0 )
2783 // makeQuadratic( occgeo.emap, _mesh );
2784 // if ( initState._nbFaces > 0 )
2785 // makeQuadratic( occgeo.fmap, _mesh );
2787 catch (Standard_Failure& ex)
2789 if ( comment.empty() ) // do not overwrite a previos error
2790 comment << "Exception in netgen at passing to 2nd order ";
2792 catch (netgen::NgException exc)
2794 if ( comment.empty() ) // do not overwrite a previos error
2795 comment << exc.What();
2800 _ticTime = 0.98 / _progressTic;
2802 int nbNod = _ngMesh->GetNP();
2803 int nbSeg = _ngMesh->GetNSeg();
2804 int nbFac = _ngMesh->GetNSE();
2805 int nbVol = _ngMesh->GetNE();
2806 bool isOK = ( !err && (_isVolume ? (nbVol > 0) : (nbFac > 0)) );
2808 MESSAGE((err ? "Mesh Generation failure" : "End of Mesh Generation") <<
2809 ", nb nodes: " << nbNod <<
2810 ", nb segments: " << nbSeg <<
2811 ", nb faces: " << nbFac <<
2812 ", nb volumes: " << nbVol);
2814 // Feed back the SMESHDS with the generated Nodes and Elements
2815 if ( true /*isOK*/ ) // get whatever built
2817 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2819 if ( quadHelper.GetIsQuadratic() ) // remove free nodes
2820 for ( size_t i = 0; i < nodeVec.size(); ++i )
2821 if ( nodeVec[i] && nodeVec[i]->NbInverseElements() == 0 )
2822 _mesh->GetMeshDS()->RemoveFreeNode( nodeVec[i], 0, /*fromGroups=*/false );
2824 SMESH_ComputeErrorPtr readErr = ReadErrors(nodeVec);
2825 if ( readErr && !readErr->myBadElements.empty() )
2828 if ( !comment.empty() && !readErr->myComment.empty() ) comment += "\n";
2829 comment += readErr->myComment;
2831 if ( error->IsOK() && ( !isOK || comment.size() > 0 ))
2832 error->myName = COMPERR_ALGO_FAILED;
2833 if ( !comment.empty() )
2834 error->myComment = comment;
2836 // SetIsAlwaysComputed( true ) to empty sub-meshes, which
2837 // appear if the geometry contains coincident sub-shape due
2838 // to bool merge_solids = 1; in netgen/libsrc/occ/occgenmesh.cpp
2839 const int nbMaps = 2;
2840 const TopTools_IndexedMapOfShape* geoMaps[nbMaps] =
2841 { & occgeo.vmap, & occgeo.emap/*, & occgeo.fmap*/ };
2842 for ( int iMap = 0; iMap < nbMaps; ++iMap )
2843 for (int i = 1; i <= geoMaps[iMap]->Extent(); i++)
2844 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( geoMaps[iMap]->FindKey(i)))
2845 if ( !sm->IsMeshComputed() )
2846 sm->SetIsAlwaysComputed( true );
2848 // set bad compute error to subshapes of all failed sub-shapes
2849 if ( !error->IsOK() )
2851 bool pb2D = false, pb3D = false;
2852 for (int i = 1; i <= occgeo.fmap.Extent(); i++) {
2853 int status = occgeo.facemeshstatus[i-1];
2854 if (status == 1 ) continue;
2855 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.fmap( i ))) {
2856 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
2857 if ( !smError || smError->IsOK() ) {
2859 smError.reset( new SMESH_ComputeError( *error ));
2861 smError.reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED, "Ignored" ));
2862 if ( SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
2863 smError->myName = COMPERR_WARNING;
2865 pb2D = pb2D || smError->IsKO();
2868 if ( !pb2D ) // all faces are OK
2869 for (int i = 1; i <= occgeo.somap.Extent(); i++)
2870 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.somap( i )))
2872 bool smComputed = nbVol && !sm->IsEmpty();
2873 if ( smComputed && internals.hasInternalVertexInSolid( sm->GetId() ))
2875 int nbIntV = internals.getSolidsWithVertices().find( sm->GetId() )->second.size();
2876 SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
2877 smComputed = ( smDS->NbElements() > 0 || smDS->NbNodes() > nbIntV );
2879 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
2880 if ( !smComputed && ( !smError || smError->IsOK() ))
2882 smError.reset( new SMESH_ComputeError( *error ));
2883 if ( nbVol && SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
2885 smError->myName = COMPERR_WARNING;
2887 else if ( !smError->myBadElements.empty() ) // bad surface mesh
2889 if ( !hasBadElemOnSolid( smError->myBadElements, sm ))
2893 pb3D = pb3D || ( smError && smError->IsKO() );
2895 if ( !pb2D && !pb3D )
2896 err = 0; // no fatal errors, only warnings
2899 ngLib._isComputeOk = !err;
2904 //=============================================================================
2908 //=============================================================================
2909 bool NETGENPlugin_Mesher::Evaluate(MapShapeNbElems& aResMap)
2911 netgen::MeshingParameters& mparams = netgen::mparam;
2914 // -------------------------
2915 // Prepare OCC geometry
2916 // -------------------------
2917 netgen::OCCGeometry occgeo;
2918 list< SMESH_subMesh* > meshedSM[4]; // for 0-3 dimensions
2919 NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
2920 PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
2922 bool tooManyElems = false;
2923 const int hugeNb = std::numeric_limits<int>::max() / 100;
2928 // pass 1D simple parameters to NETGEN
2931 // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
2932 mparams.uselocalh = false;
2933 mparams.grading = 0.8; // not limitited size growth
2935 if ( _simpleHyp->GetNumberOfSegments() )
2937 mparams.maxh = occgeo.boundingbox.Diam();
2940 mparams.maxh = _simpleHyp->GetLocalLength();
2943 if ( mparams.maxh == 0.0 )
2944 mparams.maxh = occgeo.boundingbox.Diam();
2945 if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
2946 mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
2948 // let netgen create _ngMesh and calculate element size on not meshed shapes
2949 NETGENPlugin_NetgenLibWrapper ngLib;
2950 netgen::Mesh *ngMesh = NULL;
2954 int startWith = netgen::MESHCONST_ANALYSE;
2955 int endWith = netgen::MESHCONST_MESHEDGES;
2957 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, mparams, startWith, endWith);
2959 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, startWith, endWith, optstr);
2962 if(netgen::multithread.terminate)
2965 ngLib.setMesh(( Ng_Mesh*) ngMesh );
2967 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( _shape ))
2968 sm->GetComputeError().reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED ));
2973 // Pass 1D simple parameters to NETGEN
2974 // --------------------------------
2975 int nbSeg = _simpleHyp->GetNumberOfSegments();
2976 double segSize = _simpleHyp->GetLocalLength();
2977 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
2979 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
2981 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
2982 setLocalSize( e, segSize, *ngMesh );
2985 else // if ( ! _simpleHyp )
2987 // Local size on vertices and edges
2988 // --------------------------------
2989 for(std::map<int,double>::const_iterator it=EdgeId2LocalSize.begin(); it!=EdgeId2LocalSize.end(); it++)
2991 int key = (*it).first;
2992 double hi = (*it).second;
2993 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2994 const TopoDS_Edge& e = TopoDS::Edge(shape);
2995 setLocalSize( e, hi, *ngMesh );
2997 for(std::map<int,double>::const_iterator it=VertexId2LocalSize.begin(); it!=VertexId2LocalSize.end(); it++)
2999 int key = (*it).first;
3000 double hi = (*it).second;
3001 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3002 const TopoDS_Vertex& v = TopoDS::Vertex(shape);
3003 gp_Pnt p = BRep_Tool::Pnt(v);
3004 NETGENPlugin_Mesher::RestrictLocalSize( *ngMesh, p.XYZ(), hi );
3006 for(map<int,double>::const_iterator it=FaceId2LocalSize.begin();
3007 it!=FaceId2LocalSize.end(); it++)
3009 int key = (*it).first;
3010 double val = (*it).second;
3011 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3012 int faceNgID = occgeo.fmap.FindIndex(shape);
3013 occgeo.SetFaceMaxH(faceNgID, val);
3014 for ( TopExp_Explorer edgeExp( shape, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
3015 setLocalSize( TopoDS::Edge( edgeExp.Current() ), val, *ngMesh );
3018 // calculate total nb of segments and length of edges
3019 double fullLen = 0.0;
3021 int entity = mparams.secondorder > 0 ? SMDSEntity_Quad_Edge : SMDSEntity_Edge;
3022 TopTools_DataMapOfShapeInteger Edge2NbSeg;
3023 for (TopExp_Explorer exp(_shape, TopAbs_EDGE); exp.More(); exp.Next())
3025 TopoDS_Edge E = TopoDS::Edge( exp.Current() );
3026 if( !Edge2NbSeg.Bind(E,0) )
3029 double aLen = SMESH_Algo::EdgeLength(E);
3032 vector<int>& aVec = aResMap[_mesh->GetSubMesh(E)];
3034 aVec.resize( SMDSEntity_Last, 0);
3036 fullNbSeg += aVec[ entity ];
3039 // store nb of segments computed by Netgen
3040 NCollection_Map<Link> linkMap;
3041 for (int i = 1; i <= ngMesh->GetNSeg(); ++i )
3043 const netgen::Segment& seg = ngMesh->LineSegment(i);
3044 Link link(seg[0], seg[1]);
3045 if ( !linkMap.Add( link )) continue;
3046 int aGeomEdgeInd = seg.epgeominfo[0].edgenr;
3047 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
3049 vector<int>& aVec = aResMap[_mesh->GetSubMesh(occgeo.emap(aGeomEdgeInd))];
3053 // store nb of nodes on edges computed by Netgen
3054 TopTools_DataMapIteratorOfDataMapOfShapeInteger Edge2NbSegIt(Edge2NbSeg);
3055 for (; Edge2NbSegIt.More(); Edge2NbSegIt.Next())
3057 vector<int>& aVec = aResMap[_mesh->GetSubMesh(Edge2NbSegIt.Key())];
3058 if ( aVec[ entity ] > 1 && aVec[ SMDSEntity_Node ] == 0 )
3059 aVec[SMDSEntity_Node] = mparams.secondorder > 0 ? 2*aVec[ entity ]-1 : aVec[ entity ]-1;
3061 fullNbSeg += aVec[ entity ];
3062 Edge2NbSeg( Edge2NbSegIt.Key() ) = aVec[ entity ];
3064 if ( fullNbSeg == 0 )
3071 if ( double area = _simpleHyp->GetMaxElementArea() ) {
3073 mparams.maxh = sqrt(2. * area/sqrt(3.0));
3074 mparams.grading = 0.4; // moderate size growth
3077 // length from edges
3078 mparams.maxh = fullLen/fullNbSeg;
3079 mparams.grading = 0.2; // slow size growth
3082 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3083 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3085 for (TopExp_Explorer exp(_shape, TopAbs_FACE); exp.More(); exp.Next())
3087 TopoDS_Face F = TopoDS::Face( exp.Current() );
3088 SMESH_subMesh *sm = _mesh->GetSubMesh(F);
3090 BRepGProp::SurfaceProperties(F,G);
3091 double anArea = G.Mass();
3092 tooManyElems = tooManyElems || ( anArea/hugeNb > mparams.maxh*mparams.maxh );
3094 if ( !tooManyElems )
3096 TopTools_MapOfShape egdes;
3097 for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next())
3098 if ( egdes.Add( exp1.Current() ))
3099 nb1d += Edge2NbSeg.Find(exp1.Current());
3101 int nbFaces = tooManyElems ? hugeNb : int( 4*anArea / (mparams.maxh*mparams.maxh*sqrt(3.)));
3102 int nbNodes = tooManyElems ? hugeNb : (( nbFaces*3 - (nb1d-1)*2 ) / 6 + 1 );
3104 vector<int> aVec(SMDSEntity_Last, 0);
3105 if( mparams.secondorder > 0 ) {
3106 int nb1d_in = (nbFaces*3 - nb1d) / 2;
3107 aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
3108 aVec[SMDSEntity_Quad_Triangle] = nbFaces;
3111 aVec[SMDSEntity_Node] = Max ( nbNodes, 0 );
3112 aVec[SMDSEntity_Triangle] = nbFaces;
3114 aResMap[sm].swap(aVec);
3121 // pass 3D simple parameters to NETGEN
3122 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
3123 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
3125 if ( double vol = simple3d->GetMaxElementVolume() ) {
3127 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
3128 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3131 // using previous length from faces
3133 mparams.grading = 0.4;
3134 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3137 BRepGProp::VolumeProperties(_shape,G);
3138 double aVolume = G.Mass();
3139 double tetrVol = 0.1179*mparams.maxh*mparams.maxh*mparams.maxh;
3140 tooManyElems = tooManyElems || ( aVolume/hugeNb > tetrVol );
3141 int nbVols = tooManyElems ? hugeNb : int(aVolume/tetrVol);
3142 int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
3143 vector<int> aVec(SMDSEntity_Last, 0 );
3144 if ( tooManyElems ) // avoid FPE
3146 aVec[SMDSEntity_Node] = hugeNb;
3147 aVec[ mparams.secondorder > 0 ? SMDSEntity_Quad_Tetra : SMDSEntity_Tetra] = hugeNb;
3151 if( mparams.secondorder > 0 ) {
3152 aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
3153 aVec[SMDSEntity_Quad_Tetra] = nbVols;
3156 aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
3157 aVec[SMDSEntity_Tetra] = nbVols;
3160 SMESH_subMesh *sm = _mesh->GetSubMesh(_shape);
3161 aResMap[sm].swap(aVec);
3167 double NETGENPlugin_Mesher::GetProgress(const SMESH_Algo* holder,
3168 const int * algoProgressTic,
3169 const double * algoProgress) const
3171 ((int&) _progressTic ) = *algoProgressTic + 1;
3173 if ( !_occgeom ) return 0;
3175 double progress = -1;
3178 if ( _ticTime < 0 && netgen::multithread.task[0] == 'O'/*Optimizing surface*/ )
3180 ((double&) _ticTime ) = edgeFaceMeshingTime / _totalTime / _progressTic;
3182 else if ( !_optimize /*&& _occgeom->fmap.Extent() > 1*/ )
3184 int doneShapeIndex = -1;
3185 while ( doneShapeIndex+1 < _occgeom->facemeshstatus.Size() &&
3186 _occgeom->facemeshstatus[ doneShapeIndex+1 ])
3188 if ( doneShapeIndex+1 != _curShapeIndex )
3190 ((int&) _curShapeIndex) = doneShapeIndex+1;
3191 double doneShapeRate = _curShapeIndex / double( _occgeom->fmap.Extent() );
3192 double doneTime = edgeMeshingTime + doneShapeRate * faceMeshingTime;
3193 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3194 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3195 // << " " << doneTime / _totalTime / _progressTic << endl;
3199 else if ( !_optimize && _occgeom->somap.Extent() > 1 )
3201 int curShapeIndex = _curShapeIndex;
3202 if ( _ngMesh->GetNE() > 0 )
3204 netgen::Element el = (*_ngMesh)[netgen::ElementIndex( _ngMesh->GetNE()-1 )];
3205 curShapeIndex = el.GetIndex();
3207 if ( curShapeIndex != _curShapeIndex )
3209 ((int&) _curShapeIndex) = curShapeIndex;
3210 double doneShapeRate = _curShapeIndex / double( _occgeom->somap.Extent() );
3211 double doneTime = edgeFaceMeshingTime + doneShapeRate * voluMeshingTime;
3212 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3213 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3214 // << " " << doneTime / _totalTime / _progressTic << endl;
3218 progress = Max( *algoProgressTic * _ticTime, *algoProgress );
3221 ((int&) *algoProgressTic )++;
3222 ((double&) *algoProgress) = progress;
3224 //cout << progress << " " << *algoProgressTic << " " << netgen::multithread.task << " "<< _ticTime << endl;
3226 return Min( progress, 0.99 );
3229 //================================================================================
3231 * \brief Remove "test.out" and "problemfaces" files in current directory
3233 //================================================================================
3235 void NETGENPlugin_Mesher::RemoveTmpFiles()
3237 bool rm = SMESH_File("test.out").remove() ;
3239 if (rm && netgen::testout)
3241 delete netgen::testout;
3242 netgen::testout = 0;
3245 SMESH_File("problemfaces").remove();
3246 SMESH_File("occmesh.rep").remove();
3249 //================================================================================
3251 * \brief Read mesh entities preventing successful computation from "test.out" file
3253 //================================================================================
3255 SMESH_ComputeErrorPtr
3256 NETGENPlugin_Mesher::ReadErrors(const vector<const SMDS_MeshNode* >& nodeVec)
3258 SMESH_ComputeErrorPtr err = SMESH_ComputeError::New
3259 (COMPERR_BAD_INPUT_MESH, "Some edges multiple times in surface mesh");
3260 SMESH_File file("test.out");
3262 vector<int> three1(3), three2(3);
3263 const char* badEdgeStr = " multiple times in surface mesh";
3264 const int badEdgeStrLen = strlen( badEdgeStr );
3266 while( !file.eof() )
3268 if ( strncmp( file, "Edge ", 5 ) == 0 &&
3269 file.getInts( two ) &&
3270 strncmp( file, badEdgeStr, badEdgeStrLen ) == 0 &&
3271 two[0] < nodeVec.size() && two[1] < nodeVec.size())
3273 err->myBadElements.push_back( new SMDS_LinearEdge( nodeVec[ two[0]], nodeVec[ two[1]] ));
3274 file += badEdgeStrLen;
3276 else if ( strncmp( file, "Intersecting: ", 14 ) == 0 )
3279 // openelement 18 with open element 126
3283 const char* pos = file;
3284 bool ok = ( strncmp( file, "openelement ", 12 ) == 0 );
3285 ok = ok && file.getInts( two );
3286 ok = ok && file.getInts( three1 );
3287 ok = ok && file.getInts( three2 );
3288 for ( int i = 0; ok && i < 3; ++i )
3289 ok = ( three1[i] < nodeVec.size() && nodeVec[ three1[i]]);
3290 for ( int i = 0; ok && i < 3; ++i )
3291 ok = ( three2[i] < nodeVec.size() && nodeVec[ three2[i]]);
3294 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three1[0]],
3295 nodeVec[ three1[1]],
3296 nodeVec[ three1[2]]));
3297 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three2[0]],
3298 nodeVec[ three2[1]],
3299 nodeVec[ three2[2]]));
3300 err->myComment = "Intersecting triangles";
3314 size_t nbBadElems = err->myBadElements.size();
3321 //================================================================================
3323 * \brief Write a python script creating an equivalent SALOME mesh.
3324 * This is useful to see what mesh is passed as input for the next step of mesh
3325 * generation (of mesh of higher dimension)
3327 //================================================================================
3329 void NETGENPlugin_Mesher::toPython( const netgen::Mesh* ngMesh,
3330 const std::string& pyFile)
3332 ofstream outfile(pyFile.c_str(), ios::out);
3333 if ( !outfile ) return;
3335 outfile << "import SMESH" << endl
3336 << "from salome.smesh import smeshBuilder" << endl
3337 << "smesh = smeshBuilder.New(salome.myStudy)" << endl
3338 << "mesh = smesh.Mesh()" << endl << endl;
3340 using namespace netgen;
3342 for (pi = PointIndex::BASE;
3343 pi < ngMesh->GetNP()+PointIndex::BASE; pi++)
3345 outfile << "mesh.AddNode( ";
3346 outfile << (*ngMesh)[pi](0) << ", ";
3347 outfile << (*ngMesh)[pi](1) << ", ";
3348 outfile << (*ngMesh)[pi](2) << ") ## "<< pi << endl;
3351 int nbDom = ngMesh->GetNDomains();
3352 for ( int i = 0; i < nbDom; ++i )
3353 outfile<< "grp" << i+1 << " = mesh.CreateEmptyGroup( SMESH.FACE, 'domain"<< i+1 << "')"<< endl;
3355 SurfaceElementIndex sei;
3356 for (sei = 0; sei < ngMesh->GetNSE(); sei++)
3358 outfile << "mesh.AddFace([ ";
3359 Element2d sel = (*ngMesh)[sei];
3360 for (int j = 0; j < sel.GetNP(); j++)
3361 outfile << sel[j] << ( j+1 < sel.GetNP() ? ", " : " ])");
3362 if ( sel.IsDeleted() ) outfile << " ## IsDeleted ";
3365 if ((*ngMesh)[sei].GetIndex())
3367 if ( int dom1 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainIn())
3368 outfile << "grp"<< dom1 <<".Add([ " << (int)sei+1 << " ])" << endl;
3369 if ( int dom2 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainOut())
3370 outfile << "grp"<< dom2 <<".Add([ " << (int)sei+1 << " ])" << endl;
3374 for (ElementIndex ei = 0; ei < ngMesh->GetNE(); ei++)
3376 Element el = (*ngMesh)[ei];
3377 outfile << "mesh.AddVolume([ ";
3378 for (int j = 0; j < el.GetNP(); j++)
3379 outfile << el[j] << ( j+1 < el.GetNP() ? ", " : " ])");
3383 for (int i = 1; i <= ngMesh->GetNSeg(); i++)
3385 const Segment & seg = ngMesh->LineSegment (i);
3386 outfile << "mesh.AddEdge([ "
3388 << seg[1] << " ])" << endl;
3390 cout << "Write " << pyFile << endl;
3393 //================================================================================
3395 * \brief Constructor of NETGENPlugin_ngMeshInfo
3397 //================================================================================
3399 NETGENPlugin_ngMeshInfo::NETGENPlugin_ngMeshInfo( netgen::Mesh* ngMesh):
3404 _nbNodes = ngMesh->GetNP();
3405 _nbSegments = ngMesh->GetNSeg();
3406 _nbFaces = ngMesh->GetNSE();
3407 _nbVolumes = ngMesh->GetNE();
3411 _nbNodes = _nbSegments = _nbFaces = _nbVolumes = 0;
3415 //================================================================================
3417 * \brief Copy LocalH member from one netgen mesh to another
3419 //================================================================================
3421 void NETGENPlugin_ngMeshInfo::transferLocalH( netgen::Mesh* fromMesh,
3422 netgen::Mesh* toMesh )
3424 if ( !fromMesh->LocalHFunctionGenerated() ) return;
3425 if ( !toMesh->LocalHFunctionGenerated() )
3427 toMesh->CalcLocalH(netgen::mparam.grading);
3429 toMesh->CalcLocalH();
3432 const size_t size = sizeof( netgen::LocalH );
3433 _copyOfLocalH = new char[ size ];
3434 memcpy( (void*)_copyOfLocalH, (void*)&toMesh->LocalHFunction(), size );
3435 memcpy( (void*)&toMesh->LocalHFunction(), (void*)&fromMesh->LocalHFunction(), size );
3438 //================================================================================
3440 * \brief Restore LocalH member of a netgen mesh
3442 //================================================================================
3444 void NETGENPlugin_ngMeshInfo::restoreLocalH( netgen::Mesh* toMesh )
3446 if ( _copyOfLocalH )
3448 const size_t size = sizeof( netgen::LocalH );
3449 memcpy( (void*)&toMesh->LocalHFunction(), (void*)_copyOfLocalH, size );
3450 delete [] _copyOfLocalH;
3455 //================================================================================
3457 * \brief Find "internal" sub-shapes
3459 //================================================================================
3461 NETGENPlugin_Internals::NETGENPlugin_Internals( SMESH_Mesh& mesh,
3462 const TopoDS_Shape& shape,
3464 : _mesh( mesh ), _is3D( is3D )
3466 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3468 TopExp_Explorer f,e;
3469 for ( f.Init( shape, TopAbs_FACE ); f.More(); f.Next() )
3471 int faceID = meshDS->ShapeToIndex( f.Current() );
3473 // find not computed internal edges
3475 for ( e.Init( f.Current().Oriented(TopAbs_FORWARD), TopAbs_EDGE ); e.More(); e.Next() )
3476 if ( e.Current().Orientation() == TopAbs_INTERNAL )
3478 SMESH_subMesh* eSM = mesh.GetSubMesh( e.Current() );
3479 if ( eSM->IsEmpty() )
3481 _e2face.insert( make_pair( eSM->GetId(), faceID ));
3482 for ( TopoDS_Iterator v(e.Current()); v.More(); v.Next() )
3483 _e2face.insert( make_pair( meshDS->ShapeToIndex( v.Value() ), faceID ));
3487 // find internal vertices in a face
3488 set<int> intVV; // issue 0020850 where same vertex is twice in a face
3489 for ( TopoDS_Iterator fSub( f.Current() ); fSub.More(); fSub.Next())
3490 if ( fSub.Value().ShapeType() == TopAbs_VERTEX )
3492 int vID = meshDS->ShapeToIndex( fSub.Value() );
3493 if ( intVV.insert( vID ).second )
3494 _f2v[ faceID ].push_back( vID );
3499 // find internal faces and their subshapes where nodes are to be doubled
3500 // to make a crack with non-sewed borders
3502 if ( f.Current().Orientation() == TopAbs_INTERNAL )
3504 _intShapes.insert( meshDS->ShapeToIndex( f.Current() ));
3507 list< TopoDS_Shape > edges;
3508 for ( e.Init( f.Current(), TopAbs_EDGE ); e.More(); e.Next())
3509 if ( SMESH_MesherHelper::NbAncestors( e.Current(), mesh, TopAbs_FACE ) > 1 )
3511 _intShapes.insert( meshDS->ShapeToIndex( e.Current() ));
3512 edges.push_back( e.Current() );
3513 // find border faces
3514 PShapeIteratorPtr fIt =
3515 SMESH_MesherHelper::GetAncestors( edges.back(),mesh,TopAbs_FACE );
3516 while ( const TopoDS_Shape* pFace = fIt->next() )
3517 if ( !pFace->IsSame( f.Current() ))
3518 _borderFaces.insert( meshDS->ShapeToIndex( *pFace ));
3521 // we consider vertex internal if it is shared by more than one internal edge
3522 list< TopoDS_Shape >::iterator edge = edges.begin();
3523 for ( ; edge != edges.end(); ++edge )
3524 for ( TopoDS_Iterator v( *edge ); v.More(); v.Next() )
3526 set<int> internalEdges;
3527 PShapeIteratorPtr eIt =
3528 SMESH_MesherHelper::GetAncestors( v.Value(),mesh,TopAbs_EDGE );
3529 while ( const TopoDS_Shape* pEdge = eIt->next() )
3531 int edgeID = meshDS->ShapeToIndex( *pEdge );
3532 if ( isInternalShape( edgeID ))
3533 internalEdges.insert( edgeID );
3535 if ( internalEdges.size() > 1 )
3536 _intShapes.insert( meshDS->ShapeToIndex( v.Value() ));
3540 } // loop on geom faces
3542 // find vertices internal in solids
3545 for ( TopExp_Explorer so(shape, TopAbs_SOLID); so.More(); so.Next())
3547 int soID = meshDS->ShapeToIndex( so.Current() );
3548 for ( TopoDS_Iterator soSub( so.Current() ); soSub.More(); soSub.Next())
3549 if ( soSub.Value().ShapeType() == TopAbs_VERTEX )
3550 _s2v[ soID ].push_back( meshDS->ShapeToIndex( soSub.Value() ));
3555 //================================================================================
3557 * \brief Find mesh faces on non-internal geom faces sharing internal edge
3558 * some nodes of which are to be doubled to make the second border of the "crack"
3560 //================================================================================
3562 void NETGENPlugin_Internals::findBorderElements( TIDSortedElemSet & borderElems )
3564 if ( _intShapes.empty() ) return;
3566 SMESH_Mesh& mesh = const_cast<SMESH_Mesh&>(_mesh);
3567 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3569 // loop on internal geom edges
3570 set<int>::const_iterator intShapeId = _intShapes.begin();
3571 for ( ; intShapeId != _intShapes.end(); ++intShapeId )
3573 const TopoDS_Shape& s = meshDS->IndexToShape( *intShapeId );
3574 if ( s.ShapeType() != TopAbs_EDGE ) continue;
3576 // get internal and non-internal geom faces sharing the internal edge <s>
3578 set<int>::iterator bordFace = _borderFaces.end();
3579 PShapeIteratorPtr faces = SMESH_MesherHelper::GetAncestors( s, _mesh, TopAbs_FACE );
3580 while ( const TopoDS_Shape* pFace = faces->next() )
3582 int faceID = meshDS->ShapeToIndex( *pFace );
3583 if ( isInternalShape( faceID ))
3586 bordFace = _borderFaces.insert( faceID ).first;
3588 if ( bordFace == _borderFaces.end() || !intFace ) continue;
3590 // get all links of mesh faces on internal geom face sharing nodes on edge <s>
3591 set< SMESH_OrientedLink > links; //!< links of faces on internal geom face
3592 list<const SMDS_MeshElement*> suspectFaces[2]; //!< mesh faces on border geom faces
3593 int nbSuspectFaces = 0;
3594 SMESHDS_SubMesh* intFaceSM = meshDS->MeshElements( intFace );
3595 if ( !intFaceSM || intFaceSM->NbElements() == 0 ) continue;
3596 SMESH_subMeshIteratorPtr smIt = mesh.GetSubMesh( s )->getDependsOnIterator(true,true);
3597 while ( smIt->more() )
3599 SMESHDS_SubMesh* sm = smIt->next()->GetSubMeshDS();
3600 if ( !sm ) continue;
3601 SMDS_NodeIteratorPtr nIt = sm->GetNodes();
3602 while ( nIt->more() )
3604 const SMDS_MeshNode* nOnEdge = nIt->next();
3605 SMDS_ElemIteratorPtr fIt = nOnEdge->GetInverseElementIterator(SMDSAbs_Face);
3606 while ( fIt->more() )
3608 const SMDS_MeshElement* f = fIt->next();
3609 int nbNodes = f->NbNodes() / ( f->IsQuadratic() ? 2 : 1 );
3610 if ( intFaceSM->Contains( f ))
3612 for ( int i = 0; i < nbNodes; ++i )
3613 links.insert( SMESH_OrientedLink( f->GetNode(i), f->GetNode((i+1)%nbNodes)));
3618 for ( int i = 0; i < nbNodes; ++i )
3619 nbDblNodes += isInternalShape( f->GetNode(i)->getshapeId() );
3621 suspectFaces[ nbDblNodes < 2 ].push_back( f );
3627 // suspectFaces[0] having link with same orientation as mesh faces on
3628 // the internal geom face are <borderElems>. suspectFaces[1] have
3629 // only one node on edge <s>, we decide on them later (at the 2nd loop)
3630 // by links of <borderElems> found at the 1st and 2nd loops
3631 set< SMESH_OrientedLink > borderLinks;
3632 for ( int isPostponed = 0; isPostponed < 2; ++isPostponed )
3634 list<const SMDS_MeshElement*>::iterator fIt = suspectFaces[isPostponed].begin();
3635 for ( int nbF = 0; fIt != suspectFaces[isPostponed].end(); ++fIt, ++nbF )
3637 const SMDS_MeshElement* f = *fIt;
3638 bool isBorder = false, linkFound = false, borderLinkFound = false;
3639 list< SMESH_OrientedLink > faceLinks;
3640 int nbNodes = f->NbNodes() / ( f->IsQuadratic() ? 2 : 1 );
3641 for ( int i = 0; i < nbNodes; ++i )
3643 SMESH_OrientedLink link( f->GetNode(i), f->GetNode((i+1)%nbNodes));
3644 faceLinks.push_back( link );
3647 set< SMESH_OrientedLink >::iterator foundLink = links.find( link );
3648 if ( foundLink != links.end() )
3651 isBorder = ( foundLink->_reversed == link._reversed );
3652 if ( !isBorder && !isPostponed ) break;
3653 faceLinks.pop_back();
3655 else if ( isPostponed && !borderLinkFound )
3657 foundLink = borderLinks.find( link );
3658 if ( foundLink != borderLinks.end() )
3660 borderLinkFound = true;
3661 isBorder = ( foundLink->_reversed != link._reversed );
3668 borderElems.insert( f );
3669 borderLinks.insert( faceLinks.begin(), faceLinks.end() );
3671 else if ( !linkFound && !borderLinkFound )
3673 suspectFaces[1].push_back( f );
3674 if ( nbF > 2 * nbSuspectFaces )
3675 break; // dead loop protection
3682 //================================================================================
3684 * \brief put internal shapes in maps and fill in submeshes to precompute
3686 //================================================================================
3688 void NETGENPlugin_Internals::getInternalEdges( TopTools_IndexedMapOfShape& fmap,
3689 TopTools_IndexedMapOfShape& emap,
3690 TopTools_IndexedMapOfShape& vmap,
3691 list< SMESH_subMesh* > smToPrecompute[])
3693 if ( !hasInternalEdges() ) return;
3694 map<int,int>::const_iterator ev_face = _e2face.begin();
3695 for ( ; ev_face != _e2face.end(); ++ev_face )
3697 const TopoDS_Shape& ev = _mesh.GetMeshDS()->IndexToShape( ev_face->first );
3698 const TopoDS_Shape& face = _mesh.GetMeshDS()->IndexToShape( ev_face->second );
3700 ( ev.ShapeType() == TopAbs_EDGE ? emap : vmap ).Add( ev );
3702 //cout<<"INTERNAL EDGE or VERTEX "<<ev_face->first<<" on face "<<ev_face->second<<endl;
3704 smToPrecompute[ MeshDim_1D ].push_back( _mesh.GetSubMeshContaining( ev_face->first ));
3708 //================================================================================
3710 * \brief return shapes and submeshes to be meshed and already meshed boundary submeshes
3712 //================================================================================
3714 void NETGENPlugin_Internals::getInternalFaces( TopTools_IndexedMapOfShape& fmap,
3715 TopTools_IndexedMapOfShape& emap,
3716 list< SMESH_subMesh* >& intFaceSM,
3717 list< SMESH_subMesh* >& boundarySM)
3719 if ( !hasInternalFaces() ) return;
3721 // <fmap> and <emap> are for not yet meshed shapes
3722 // <intFaceSM> is for submeshes of faces
3723 // <boundarySM> is for meshed edges and vertices
3728 set<int> shapeIDs ( _intShapes );
3729 if ( !_borderFaces.empty() )
3730 shapeIDs.insert( _borderFaces.begin(), _borderFaces.end() );
3732 set<int>::const_iterator intS = shapeIDs.begin();
3733 for ( ; intS != shapeIDs.end(); ++intS )
3735 SMESH_subMesh* sm = _mesh.GetSubMeshContaining( *intS );
3737 if ( sm->GetSubShape().ShapeType() != TopAbs_FACE ) continue;
3739 intFaceSM.push_back( sm );
3741 // add submeshes of not computed internal faces
3742 if ( !sm->IsEmpty() ) continue;
3744 SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(true,true);
3745 while ( smIt->more() )
3748 const TopoDS_Shape& s = sm->GetSubShape();
3750 if ( sm->IsEmpty() )
3753 switch ( s.ShapeType() ) {
3754 case TopAbs_FACE: fmap.Add ( s ); break;
3755 case TopAbs_EDGE: emap.Add ( s ); break;
3761 if ( s.ShapeType() != TopAbs_FACE )
3762 boundarySM.push_back( sm );
3768 //================================================================================
3770 * \brief Return true if given shape is to be precomputed in order to be correctly
3771 * added to netgen mesh
3773 //================================================================================
3775 bool NETGENPlugin_Internals::isShapeToPrecompute(const TopoDS_Shape& s)
3777 int shapeID = _mesh.GetMeshDS()->ShapeToIndex( s );
3778 switch ( s.ShapeType() ) {
3779 case TopAbs_FACE : break; //return isInternalShape( shapeID ) || isBorderFace( shapeID );
3780 case TopAbs_EDGE : return isInternalEdge( shapeID );
3781 case TopAbs_VERTEX: break;
3787 //================================================================================
3789 * \brief Return SMESH
3791 //================================================================================
3793 SMESH_Mesh& NETGENPlugin_Internals::getMesh() const
3795 return const_cast<SMESH_Mesh&>( _mesh );
3798 //================================================================================
3800 * \brief Initialize netgen library
3802 //================================================================================
3804 NETGENPlugin_NetgenLibWrapper::NETGENPlugin_NetgenLibWrapper()
3808 _isComputeOk = false;
3810 if ( !getenv( "KEEP_NETGEN_OUTPUT" ))
3812 // redirect all netgen output (mycout,myerr,cout) to _outputFileName
3813 _outputFileName = getOutputFileName();
3814 netgen::mycout = new ofstream ( _outputFileName.c_str() );
3815 netgen::myerr = netgen::mycout;
3816 _coutBuffer = std::cout.rdbuf();
3818 cout << "NOTE: netgen output is redirected to file " << _outputFileName << endl;
3820 std::cout.rdbuf( netgen::mycout->rdbuf() );
3824 _ngMesh = Ng_NewMesh();
3827 //================================================================================
3829 * \brief Finish using netgen library
3831 //================================================================================
3833 NETGENPlugin_NetgenLibWrapper::~NETGENPlugin_NetgenLibWrapper()
3835 Ng_DeleteMesh( _ngMesh );
3837 NETGENPlugin_Mesher::RemoveTmpFiles();
3839 std::cout.rdbuf( _coutBuffer );
3846 //================================================================================
3848 * \brief Set netgen mesh to delete at destruction
3850 //================================================================================
3852 void NETGENPlugin_NetgenLibWrapper::setMesh( Ng_Mesh* mesh )
3855 Ng_DeleteMesh( _ngMesh );
3859 //================================================================================
3861 * \brief Return a unique file name
3863 //================================================================================
3865 std::string NETGENPlugin_NetgenLibWrapper::getOutputFileName()
3867 std::string aTmpDir = SALOMEDS_Tool::GetTmpDir();
3869 TCollection_AsciiString aGenericName = (char*)aTmpDir.c_str();
3870 aGenericName += "NETGEN_";
3872 aGenericName += getpid();
3874 aGenericName += _getpid();
3876 aGenericName += "_";
3877 aGenericName += Abs((Standard_Integer)(long) aGenericName.ToCString());
3878 aGenericName += ".out";
3880 return aGenericName.ToCString();
3883 //================================================================================
3885 * \brief Remove file with netgen output
3887 //================================================================================
3889 void NETGENPlugin_NetgenLibWrapper::removeOutputFile()
3891 if ( !_outputFileName.empty() )
3893 if ( netgen::mycout )
3895 delete netgen::mycout;
3899 string tmpDir = SALOMEDS_Tool::GetDirFromPath ( _outputFileName );
3900 string aFileName = SALOMEDS_Tool::GetNameFromPath( _outputFileName ) + ".out";
3901 SALOMEDS::ListOfFileNames_var aFiles = new SALOMEDS::ListOfFileNames;
3903 aFiles[0] = aFileName.c_str();
3905 SALOMEDS_Tool::RemoveTemporaryFiles( tmpDir.c_str(), aFiles.in(), true );