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 //=============================================================================
215 void SetLocalSize(TopoDS_Shape GeomShape, double LocalSize)
217 TopAbs_ShapeEnum GeomType = GeomShape.ShapeType();
218 if (GeomType == TopAbs_COMPOUND) {
219 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()) {
220 SetLocalSize(it.Value(), LocalSize);
225 if (! ShapesWithLocalSize.Contains(GeomShape))
226 key = ShapesWithLocalSize.Add(GeomShape);
228 key = ShapesWithLocalSize.FindIndex(GeomShape);
229 if (GeomType == TopAbs_VERTEX) {
230 VertexId2LocalSize[key] = LocalSize;
231 } else if (GeomType == TopAbs_EDGE) {
232 EdgeId2LocalSize[key] = LocalSize;
233 } else if (GeomType == TopAbs_FACE) {
234 FaceId2LocalSize[key] = LocalSize;
238 //=============================================================================
240 * Pass parameters to NETGEN
242 //=============================================================================
243 void NETGENPlugin_Mesher::SetParameters(const NETGENPlugin_Hypothesis* hyp)
247 netgen::MeshingParameters& mparams = netgen::mparam;
248 // Initialize global NETGEN parameters:
249 // maximal mesh segment size
250 mparams.maxh = hyp->GetMaxSize();
251 // maximal mesh element linear size
252 mparams.minh = hyp->GetMinSize();
253 // minimal number of segments per edge
254 mparams.segmentsperedge = hyp->GetNbSegPerEdge();
255 // rate of growth of size between elements
256 mparams.grading = hyp->GetGrowthRate();
257 // safety factor for curvatures (elements per radius)
258 mparams.curvaturesafety = hyp->GetNbSegPerRadius();
259 // create elements of second order
260 mparams.secondorder = hyp->GetSecondOrder() ? 1 : 0;
261 // quad-dominated surface meshing
262 // only triangles are allowed for volumic mesh (before realizing IMP 0021676)
264 mparams.quad = hyp->GetQuadAllowed() ? 1 : 0;
265 _optimize = hyp->GetOptimize();
266 _fineness = hyp->GetFineness();
267 mparams.uselocalh = hyp->GetSurfaceCurvature();
268 netgen::merge_solids = hyp->GetFuseEdges();
271 SMESH_Gen_i* smeshGen_i = SMESH_Gen_i::GetSMESHGen();
272 CORBA::Object_var anObject = smeshGen_i->GetNS()->Resolve("/myStudyManager");
273 SALOMEDS::StudyManager_var aStudyMgr = SALOMEDS::StudyManager::_narrow(anObject);
274 SALOMEDS::Study_var myStudy = aStudyMgr->GetStudyByID(hyp->GetStudyId());
276 const NETGENPlugin_Hypothesis::TLocalSize localSizes = hyp->GetLocalSizesAndEntries();
277 NETGENPlugin_Hypothesis::TLocalSize::const_iterator it = localSizes.begin();
278 for (it ; it != localSizes.end() ; it++)
280 std::string entry = (*it).first;
281 double val = (*it).second;
283 GEOM::GEOM_Object_var aGeomObj;
284 TopoDS_Shape S = TopoDS_Shape();
285 SALOMEDS::SObject_var aSObj = myStudy->FindObjectID( entry.c_str() );
286 if (!aSObj->_is_nil()) {
287 CORBA::Object_var obj = aSObj->GetObject();
288 aGeomObj = GEOM::GEOM_Object::_narrow(obj);
291 if ( !aGeomObj->_is_nil() )
292 S = smeshGen_i->GeomObjectToShape( aGeomObj.in() );
294 SetLocalSize(S, val);
299 //=============================================================================
301 * Pass simple parameters to NETGEN
303 //=============================================================================
305 void NETGENPlugin_Mesher::SetParameters(const NETGENPlugin_SimpleHypothesis_2D* hyp)
309 SetDefaultParameters();
312 //=============================================================================
314 * Link - a pair of integer numbers
316 //=============================================================================
320 Link(int _n1, int _n2) : n1(_n1), n2(_n2) {}
321 Link() : n1(0), n2(0) {}
324 int HashCode(const Link& aLink, int aLimit)
326 return HashCode(aLink.n1 + aLink.n2, aLimit);
329 Standard_Boolean IsEqual(const Link& aLink1, const Link& aLink2)
331 return (aLink1.n1 == aLink2.n1 && aLink1.n2 == aLink2.n2 ||
332 aLink1.n1 == aLink2.n2 && aLink1.n2 == aLink2.n1);
337 //================================================================================
339 * \brief return id of netgen point corresponding to SMDS node
341 //================================================================================
342 typedef map< const SMDS_MeshNode*, int > TNode2IdMap;
344 int ngNodeId( const SMDS_MeshNode* node,
345 netgen::Mesh& ngMesh,
346 TNode2IdMap& nodeNgIdMap)
348 int newNgId = ngMesh.GetNP() + 1;
350 TNode2IdMap::iterator node_id = nodeNgIdMap.insert( make_pair( node, newNgId )).first;
352 if ( node_id->second == newNgId)
354 #if defined(DUMP_SEGMENTS) || defined(DUMP_TRIANGLES)
355 cout << "Ng " << newNgId << " - " << node;
357 netgen::MeshPoint p( netgen::Point<3> (node->X(), node->Y(), node->Z()) );
358 ngMesh.AddPoint( p );
360 return node_id->second;
363 //================================================================================
365 * \brief Return computed EDGEs connected to the given one
367 //================================================================================
369 list< TopoDS_Edge > getConnectedEdges( const TopoDS_Edge& edge,
370 const TopoDS_Face& face,
371 const set< SMESH_subMesh* > & computedSM,
372 const SMESH_MesherHelper& helper,
373 map< SMESH_subMesh*, set< int > >& addedEdgeSM2Faces)
376 list< TopoDS_Edge > edges;
377 list< int > nbEdgesInWire;
378 int nbWires = SMESH_Block::GetOrderedEdges( face, edges, nbEdgesInWire);
380 // find <edge> within <edges>
381 list< TopoDS_Edge >::iterator eItFwd = edges.begin();
382 for ( ; eItFwd != edges.end(); ++eItFwd )
383 if ( edge.IsSame( *eItFwd ))
385 if ( eItFwd == edges.end()) return list< TopoDS_Edge>();
387 if ( eItFwd->Orientation() >= TopAbs_INTERNAL )
389 // connected INTERNAL edges returned from GetOrderedEdges() are wrongly oriented
390 // so treat each INTERNAL edge separately
391 TopoDS_Edge e = *eItFwd;
393 edges.push_back( e );
397 // get all computed EDGEs connected to <edge>
399 list< TopoDS_Edge >::iterator eItBack = eItFwd, ePrev;
400 TopoDS_Vertex vCommon;
401 TopTools_MapOfShape eAdded; // map used not to add a seam edge twice to <edges>
404 // put edges before <edge> to <edges> back
405 while ( edges.begin() != eItFwd )
406 edges.splice( edges.end(), edges, edges.begin() );
410 while ( ++eItFwd != edges.end() )
412 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( *eItFwd );
414 bool connected = TopExp::CommonVertex( *ePrev, *eItFwd, vCommon );
415 bool computed = sm->IsMeshComputed();
416 bool added = addedEdgeSM2Faces[ sm ].count( helper.GetSubShapeID() );
417 bool doubled = !eAdded.Add( *eItFwd );
418 bool orientOK = (( ePrev ->Orientation() < TopAbs_INTERNAL ) ==
419 ( eItFwd->Orientation() < TopAbs_INTERNAL ) );
420 if ( !connected || !computed || !orientOK || added || doubled )
422 // stop advancement; move edges from tail to head
423 while ( edges.back() != *ePrev )
424 edges.splice( edges.begin(), edges, --edges.end() );
430 while ( eItBack != edges.begin() )
434 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( *eItBack );
436 bool connected = TopExp::CommonVertex( *ePrev, *eItBack, vCommon );
437 bool computed = sm->IsMeshComputed();
438 bool added = addedEdgeSM2Faces[ sm ].count( helper.GetSubShapeID() );
439 bool doubled = !eAdded.Add( *eItBack );
440 bool orientOK = (( ePrev ->Orientation() < TopAbs_INTERNAL ) ==
441 ( eItBack->Orientation() < TopAbs_INTERNAL ) );
442 if ( !connected || !computed || !orientOK || added || doubled)
445 edges.erase( edges.begin(), ePrev );
449 if ( edges.front() != edges.back() )
451 // assure that the 1st vertex is meshed
452 TopoDS_Edge eLast = edges.back();
453 while ( !SMESH_Algo::VertexNode( SMESH_MesherHelper::IthVertex( 0, edges.front()), helper.GetMeshDS())
455 edges.front() != eLast )
456 edges.splice( edges.end(), edges, edges.begin() );
461 //================================================================================
463 * \brief Make triangulation of a shape precise enough
465 //================================================================================
467 void updateTriangulation( const TopoDS_Shape& shape )
469 // static set< Poly_Triangulation* > updated;
471 // TopLoc_Location loc;
472 // TopExp_Explorer fExp( shape, TopAbs_FACE );
473 // for ( ; fExp.More(); fExp.Next() )
475 // Handle(Poly_Triangulation) triangulation =
476 // BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
477 // if ( triangulation.IsNull() ||
478 // updated.insert( triangulation.operator->() ).second )
480 // BRepTools::Clean (shape);
483 BRepMesh_IncrementalMesh e(shape, 0.01, true);
485 catch (Standard_Failure)
488 // updated.erase( triangulation.operator->() );
489 // triangulation = BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
490 // updated.insert( triangulation.operator->() );
494 //================================================================================
496 * \brief Returns a medium node either existing in SMESH of created by NETGEN
497 * \param [in] corner1 - corner node 1
498 * \param [in] corner2 - corner node 2
499 * \param [in] defaultMedium - the node created by NETGEN
500 * \param [in] helper - holder of medium nodes existing in SMESH
501 * \return const SMDS_MeshNode* - the result node
503 //================================================================================
505 const SMDS_MeshNode* mediumNode( const SMDS_MeshNode* corner1,
506 const SMDS_MeshNode* corner2,
507 const SMDS_MeshNode* defaultMedium,
508 const SMESH_MesherHelper* helper)
512 TLinkNodeMap::const_iterator l2n =
513 helper->GetTLinkNodeMap().find( SMESH_TLink( corner1, corner2 ));
514 if ( l2n != helper->GetTLinkNodeMap().end() )
515 defaultMedium = l2n->second;
517 return defaultMedium;
520 //================================================================================
522 * \brief Assure that mesh on given shapes is quadratic
524 //================================================================================
526 void makeQuadratic( const TopTools_IndexedMapOfShape& shapes,
529 for ( int i = 1; i <= shapes.Extent(); ++i )
531 SMESHDS_SubMesh* smDS = mesh->GetMeshDS()->MeshElements( shapes(i) );
532 if ( !smDS ) continue;
533 SMDS_ElemIteratorPtr elemIt = smDS->GetElements();
534 if ( !elemIt->more() ) continue;
535 const SMDS_MeshElement* e = elemIt->next();
536 if ( !e || e->IsQuadratic() )
539 TIDSortedElemSet elems;
541 while ( elemIt->more() )
542 elems.insert( elems.end(), elemIt->next() );
544 SMESH_MeshEditor( mesh ).ConvertToQuadratic( /*3d=*/false, elems, /*biQuad=*/false );
550 //================================================================================
552 * \brief Initialize netgen::OCCGeometry with OCCT shape
554 //================================================================================
556 void NETGENPlugin_Mesher::PrepareOCCgeometry(netgen::OCCGeometry& occgeo,
557 const TopoDS_Shape& shape,
559 list< SMESH_subMesh* > * meshedSM,
560 NETGENPlugin_Internals* intern)
562 updateTriangulation( shape );
565 BRepBndLib::Add (shape, bb);
566 double x1,y1,z1,x2,y2,z2;
567 bb.Get (x1,y1,z1,x2,y2,z2);
568 MESSAGE("shape bounding box:\n" <<
569 "(" << x1 << " " << y1 << " " << z1 << ") " <<
570 "(" << x2 << " " << y2 << " " << z2 << ")");
571 netgen::Point<3> p1 = netgen::Point<3> (x1,y1,z1);
572 netgen::Point<3> p2 = netgen::Point<3> (x2,y2,z2);
573 occgeo.boundingbox = netgen::Box<3> (p1,p2);
575 occgeo.shape = shape;
578 // fill maps of shapes of occgeo with not yet meshed subshapes
580 // get root submeshes
581 list< SMESH_subMesh* > rootSM;
582 const int shapeID = mesh.GetMeshDS()->ShapeToIndex( shape );
583 if ( shapeID > 0 ) { // SMESH_subMesh with ID 0 may exist, don't use it!
584 rootSM.push_back( mesh.GetSubMesh( shape ));
587 for ( TopoDS_Iterator it( shape ); it.More(); it.Next() )
588 rootSM.push_back( mesh.GetSubMesh( it.Value() ));
591 // add subshapes of empty submeshes
592 list< SMESH_subMesh* >::iterator rootIt = rootSM.begin(), rootEnd = rootSM.end();
593 for ( ; rootIt != rootEnd; ++rootIt ) {
594 SMESH_subMesh * root = *rootIt;
595 SMESH_subMeshIteratorPtr smIt = root->getDependsOnIterator(/*includeSelf=*/true,
596 /*complexShapeFirst=*/true);
597 // to find a right orientation of subshapes (PAL20462)
598 TopTools_IndexedMapOfShape subShapes;
599 TopExp::MapShapes(root->GetSubShape(), subShapes);
600 while ( smIt->more() )
602 SMESH_subMesh* sm = smIt->next();
603 TopoDS_Shape shape = sm->GetSubShape();
604 if ( intern && intern->isShapeToPrecompute( shape ))
606 if ( !meshedSM || sm->IsEmpty() )
608 if ( shape.ShapeType() != TopAbs_VERTEX )
609 shape = subShapes( subShapes.FindIndex( shape ));// shape -> index -> oriented shape
610 if ( shape.Orientation() >= TopAbs_INTERNAL )
611 shape.Orientation( TopAbs_FORWARD ); // isuue 0020676
612 switch ( shape.ShapeType() ) {
613 case TopAbs_FACE : occgeo.fmap.Add( shape ); break;
614 case TopAbs_EDGE : occgeo.emap.Add( shape ); break;
615 case TopAbs_VERTEX: occgeo.vmap.Add( shape ); break;
616 case TopAbs_SOLID :occgeo.somap.Add( shape ); break;
620 // collect submeshes of meshed shapes
623 const int dim = SMESH_Gen::GetShapeDim( shape );
624 meshedSM[ dim ].push_back( sm );
628 occgeo.facemeshstatus.SetSize (occgeo.fmap.Extent());
629 occgeo.facemeshstatus = 0;
630 occgeo.face_maxh_modified.SetSize(occgeo.fmap.Extent());
631 occgeo.face_maxh_modified = 0;
632 occgeo.face_maxh.SetSize(occgeo.fmap.Extent());
633 occgeo.face_maxh = netgen::mparam.maxh;
636 //================================================================================
638 * \brief Return a default min size value suitable for the given geometry.
640 //================================================================================
642 double NETGENPlugin_Mesher::GetDefaultMinSize(const TopoDS_Shape& geom,
643 const double maxSize)
645 updateTriangulation( geom );
649 const int* pi[4] = { &i1, &i2, &i3, &i1 };
652 TopExp_Explorer fExp( geom, TopAbs_FACE );
653 for ( ; fExp.More(); fExp.Next() )
655 Handle(Poly_Triangulation) triangulation =
656 BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
657 if ( triangulation.IsNull() ) continue;
658 const double fTol = BRep_Tool::Tolerance( TopoDS::Face( fExp.Current() ));
659 const TColgp_Array1OfPnt& points = triangulation->Nodes();
660 const Poly_Array1OfTriangle& trias = triangulation->Triangles();
661 for ( int iT = trias.Lower(); iT <= trias.Upper(); ++iT )
663 trias(iT).Get( i1, i2, i3 );
664 for ( int j = 0; j < 3; ++j )
666 double dist2 = points(*pi[j]).SquareDistance( points( *pi[j+1] ));
667 if ( dist2 < minh && fTol*fTol < dist2 )
669 bb.Add( points(*pi[j]));
673 if ( minh > 0.25 * bb.SquareExtent() ) // simple geometry, rough triangulation
675 minh = 1e-3 * sqrt( bb.SquareExtent());
676 //cout << "BND BOX minh = " <<minh << endl;
680 minh = 3 * sqrt( minh ); // triangulation for visualization is rather fine
681 //cout << "TRIANGULATION minh = " <<minh << endl;
683 if ( minh > 0.5 * maxSize )
689 //================================================================================
691 * \brief Restrict size of elements at a given point
693 //================================================================================
695 void NETGENPlugin_Mesher::RestrictLocalSize(netgen::Mesh& ngMesh,
698 const bool overrideMinH)
700 if ( overrideMinH && netgen::mparam.minh > size )
702 ngMesh.SetMinimalH( size );
703 netgen::mparam.minh = size;
705 netgen::Point3d pi(p.X(), p.Y(), p.Z());
706 ngMesh.RestrictLocalH( pi, size );
709 //================================================================================
711 * \brief fill ngMesh with nodes and elements of computed submeshes
713 //================================================================================
715 bool NETGENPlugin_Mesher::FillNgMesh(netgen::OCCGeometry& occgeom,
716 netgen::Mesh& ngMesh,
717 vector<const SMDS_MeshNode*>& nodeVec,
718 const list< SMESH_subMesh* > & meshedSM,
719 SMESH_MesherHelper* quadHelper,
720 SMESH_ProxyMesh::Ptr proxyMesh)
722 TNode2IdMap nodeNgIdMap;
723 for ( int i = 1; i < nodeVec.size(); ++i )
724 nodeNgIdMap.insert( make_pair( nodeVec[i], i ));
726 TopTools_MapOfShape visitedShapes;
727 map< SMESH_subMesh*, set< int > > visitedEdgeSM2Faces;
728 set< SMESH_subMesh* > computedSM( meshedSM.begin(), meshedSM.end() );
730 SMESH_MesherHelper helper (*_mesh);
732 int faceNgID = ngMesh.GetNFD();
734 list< SMESH_subMesh* >::const_iterator smIt, smEnd = meshedSM.end();
735 for ( smIt = meshedSM.begin(); smIt != smEnd; ++smIt )
737 SMESH_subMesh* sm = *smIt;
738 if ( !visitedShapes.Add( sm->GetSubShape() ))
741 const SMESHDS_SubMesh * smDS = sm->GetSubMeshDS();
742 if ( !smDS ) continue;
744 switch ( sm->GetSubShape().ShapeType() )
746 case TopAbs_EDGE: { // EDGE
747 // ----------------------
748 TopoDS_Edge geomEdge = TopoDS::Edge( sm->GetSubShape() );
749 if ( geomEdge.Orientation() >= TopAbs_INTERNAL )
750 geomEdge.Orientation( TopAbs_FORWARD ); // issue 0020676
752 // Add ng segments for each not meshed FACE the EDGE bounds
753 PShapeIteratorPtr fIt = helper.GetAncestors( geomEdge, *sm->GetFather(), TopAbs_FACE );
754 while ( const TopoDS_Shape * anc = fIt->next() )
756 faceNgID = occgeom.fmap.FindIndex( *anc );
758 continue; // meshed face
760 int faceSMDSId = helper.GetMeshDS()->ShapeToIndex( *anc );
761 if ( visitedEdgeSM2Faces[ sm ].count( faceSMDSId ))
762 continue; // already treated EDGE
764 TopoDS_Face face = TopoDS::Face( occgeom.fmap( faceNgID ));
765 if ( face.Orientation() >= TopAbs_INTERNAL )
766 face.Orientation( TopAbs_FORWARD ); // issue 0020676
768 // get all meshed EDGEs of the FACE connected to geomEdge (issue 0021140)
769 helper.SetSubShape( face );
770 list< TopoDS_Edge > edges = getConnectedEdges( geomEdge, face, computedSM, helper,
771 visitedEdgeSM2Faces );
773 continue; // wrong ancestor?
775 // find out orientation of <edges> within <face>
776 TopoDS_Edge eNotSeam = edges.front();
777 if ( helper.HasSeam() )
779 list< TopoDS_Edge >::iterator eIt = edges.begin();
780 while ( helper.IsRealSeam( *eIt )) ++eIt;
781 if ( eIt != edges.end() )
784 TopAbs_Orientation fOri = helper.GetSubShapeOri( face, eNotSeam );
785 bool isForwad = ( fOri == eNotSeam.Orientation() || fOri >= TopAbs_INTERNAL );
787 // get all nodes from connected <edges>
788 bool isQuad = smDS->NbElements() ? smDS->GetElements()->next()->IsQuadratic() : false;
789 StdMeshers_FaceSide fSide( face, edges, _mesh, isForwad, isQuad );
790 const vector<UVPtStruct>& points = fSide.GetUVPtStruct();
791 int i, nbSeg = fSide.NbSegments();
793 // remember EDGEs of fSide to treat only once
794 for ( int iE = 0; iE < fSide.NbEdges(); ++iE )
795 visitedEdgeSM2Faces[ helper.GetMesh()->GetSubMesh( fSide.Edge(iE )) ].insert(faceSMDSId);
797 double otherSeamParam = 0;
802 int prevNgId = ngNodeId( points[0].node, ngMesh, nodeNgIdMap );
804 for ( i = 0; i < nbSeg; ++i )
806 const UVPtStruct& p1 = points[ i ];
807 const UVPtStruct& p2 = points[ i+1 ];
809 if ( p1.node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX ) //an EDGE begins
812 if ( helper.IsRealSeam( p1.node->getshapeId() ))
814 TopoDS_Edge e = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
815 isSeam = helper.IsRealSeam( e );
818 otherSeamParam = helper.GetOtherParam( helper.GetPeriodicIndex() & 1 ? p2.u : p2.v );
825 seg[1] = prevNgId = ngNodeId( p2.node, ngMesh, nodeNgIdMap );
826 // node param on curve
827 seg.epgeominfo[ 0 ].dist = p1.param;
828 seg.epgeominfo[ 1 ].dist = p2.param;
830 seg.epgeominfo[ 0 ].u = p1.u;
831 seg.epgeominfo[ 0 ].v = p1.v;
832 seg.epgeominfo[ 1 ].u = p2.u;
833 seg.epgeominfo[ 1 ].v = p2.v;
835 //geomEdge = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
836 //seg.epgeominfo[ 0 ].edgenr = seg.epgeominfo[ 1 ].edgenr = occgeom.emap.FindIndex( geomEdge );
838 //seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
839 seg.si = faceNgID; // = geom.fmap.FindIndex (face);
840 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
841 ngMesh.AddSegment (seg);
843 SMESH_TNodeXYZ np1( p1.node ), np2( p2.node );
844 RestrictLocalSize( ngMesh, 0.5*(np1+np2), (np1-np2).Modulus() );
847 cout << "Segment: " << seg.edgenr << " on SMESH face " << helper.GetMeshDS()->ShapeToIndex( face ) << endl
848 << "\tface index: " << seg.si << endl
849 << "\tp1: " << seg[0] << endl
850 << "\tp2: " << seg[1] << endl
851 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
852 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
853 //<< "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
854 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
855 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl;
856 //<< "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
860 if ( helper.GetPeriodicIndex() && 1 ) {
861 seg.epgeominfo[ 0 ].u = otherSeamParam;
862 seg.epgeominfo[ 1 ].u = otherSeamParam;
863 swap (seg.epgeominfo[0].v, seg.epgeominfo[1].v);
865 seg.epgeominfo[ 0 ].v = otherSeamParam;
866 seg.epgeominfo[ 1 ].v = otherSeamParam;
867 swap (seg.epgeominfo[0].u, seg.epgeominfo[1].u);
869 swap (seg[0], seg[1]);
870 swap (seg.epgeominfo[0].dist, seg.epgeominfo[1].dist);
871 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
872 ngMesh.AddSegment (seg);
874 cout << "Segment: " << seg.edgenr << endl
875 << "\t is SEAM (reverse) of the previous. "
876 << " Other " << (helper.GetPeriodicIndex() && 1 ? "U" : "V")
877 << " = " << otherSeamParam << endl;
880 else if ( fOri == TopAbs_INTERNAL )
882 swap (seg[0], seg[1]);
883 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
884 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
885 ngMesh.AddSegment (seg);
887 cout << "Segment: " << seg.edgenr << endl << "\t is REVERSE of the previous" << endl;
891 } // loop on geomEdge ancestors
893 if ( quadHelper ) // remember medium nodes of sub-meshes
895 SMDS_ElemIteratorPtr edges = smDS->GetElements();
896 while ( edges->more() )
898 const SMDS_MeshElement* e = edges->next();
899 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshEdge*>( e )))
905 } // case TopAbs_EDGE
907 case TopAbs_FACE: { // FACE
908 // ----------------------
909 const TopoDS_Face& geomFace = TopoDS::Face( sm->GetSubShape() );
910 helper.SetSubShape( geomFace );
911 bool isInternalFace = ( geomFace.Orientation() == TopAbs_INTERNAL );
913 // Find solids the geomFace bounds
914 int solidID1 = 0, solidID2 = 0;
915 StdMeshers_QuadToTriaAdaptor* quadAdaptor =
916 dynamic_cast<StdMeshers_QuadToTriaAdaptor*>( proxyMesh.get() );
919 solidID1 = occgeom.somap.FindIndex( quadAdaptor->GetShape() );
923 PShapeIteratorPtr solidIt = helper.GetAncestors( geomFace, *sm->GetFather(), TopAbs_SOLID);
924 while ( const TopoDS_Shape * solid = solidIt->next() )
926 int id = occgeom.somap.FindIndex ( *solid );
927 if ( solidID1 && id != solidID1 ) solidID2 = id;
931 // Add ng face descriptors of meshed faces
933 ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(faceNgID, solidID1, solidID2, 0));
935 // if second oreder is required, even already meshed faces must be passed to NETGEN
936 int fID = occgeom.fmap.Add( geomFace );
937 while ( fID < faceNgID ) // geomFace is already in occgeom.fmap, add a copy
938 fID = occgeom.fmap.Add( BRepBuilderAPI_Copy( geomFace, /*copyGeom=*/false ));
939 // Problem with the second order in a quadrangular mesh remains.
940 // 1) All quadrangles generated by NETGEN are moved to an inexistent face
941 // by FillSMesh() (find "AddFaceDescriptor")
942 // 2) Temporary triangles generated by StdMeshers_QuadToTriaAdaptor
943 // are on faces where quadrangles were.
944 // Due to these 2 points, wrong geom faces are used while conversion to qudratic
945 // of the mentioned above quadrangles and triangles
947 // Orient the face correctly in solidID1 (issue 0020206)
948 bool reverse = false;
950 TopoDS_Shape solid = occgeom.somap( solidID1 );
951 TopAbs_Orientation faceOriInSolid = helper.GetSubShapeOri( solid, geomFace );
952 if ( faceOriInSolid >= 0 )
954 helper.IsReversedSubMesh( TopoDS::Face( geomFace.Oriented( faceOriInSolid )));
957 // Add surface elements
959 netgen::Element2d tri(3);
960 tri.SetIndex ( faceNgID );
963 #ifdef DUMP_TRIANGLES
964 cout << "SMESH face " << helper.GetMeshDS()->ShapeToIndex( geomFace )
965 << " internal="<<isInternalFace << endl;
968 smDS = proxyMesh->GetSubMesh( geomFace );
970 SMDS_ElemIteratorPtr faces = smDS->GetElements();
971 while ( faces->more() )
973 const SMDS_MeshElement* f = faces->next();
974 if ( f->NbNodes() % 3 != 0 ) // not triangle
976 PShapeIteratorPtr solidIt=helper.GetAncestors(geomFace,*sm->GetFather(),TopAbs_SOLID);
977 if ( const TopoDS_Shape * solid = solidIt->next() )
978 sm = _mesh->GetSubMesh( *solid );
979 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
980 smError.reset( new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH,"Not triangle submesh"));
981 smError->myBadElements.push_back( f );
985 for ( int i = 0; i < 3; ++i )
987 const SMDS_MeshNode* node = f->GetNode( i ), * inFaceNode=0;
989 // get node UV on face
990 int shapeID = node->getshapeId();
991 if ( helper.IsSeamShape( shapeID ))
992 if ( helper.IsSeamShape( f->GetNodeWrap( i+1 )->getshapeId() ))
993 inFaceNode = f->GetNodeWrap( i-1 );
995 inFaceNode = f->GetNodeWrap( i+1 );
996 gp_XY uv = helper.GetNodeUV( geomFace, node, inFaceNode );
998 int ind = reverse ? 3-i : i+1;
999 tri.GeomInfoPi(ind).u = uv.X();
1000 tri.GeomInfoPi(ind).v = uv.Y();
1001 tri.PNum (ind) = ngNodeId( node, ngMesh, nodeNgIdMap );
1004 ngMesh.AddSurfaceElement (tri);
1005 #ifdef DUMP_TRIANGLES
1006 cout << tri << endl;
1009 if ( isInternalFace )
1011 swap( tri[1], tri[2] );
1012 ngMesh.AddSurfaceElement (tri);
1013 #ifdef DUMP_TRIANGLES
1014 cout << tri << endl;
1019 if ( quadHelper ) // remember medium nodes of sub-meshes
1021 SMDS_ElemIteratorPtr faces = smDS->GetElements();
1022 while ( faces->more() )
1024 const SMDS_MeshElement* f = faces->next();
1025 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshFace*>( f )))
1031 } // case TopAbs_FACE
1033 case TopAbs_VERTEX: { // VERTEX
1034 // --------------------------
1035 // issue 0021405. Add node only if a VERTEX is shared by a not meshed EDGE,
1036 // else netgen removes a free node and nodeVector becomes invalid
1037 PShapeIteratorPtr ansIt = helper.GetAncestors( sm->GetSubShape(),
1041 while ( const TopoDS_Shape* e = ansIt->next() )
1043 SMESH_subMesh* eSub = helper.GetMesh()->GetSubMesh( *e );
1044 if (( toAdd = ( eSub->IsEmpty() && !SMESH_Algo::isDegenerated( TopoDS::Edge( *e )))))
1049 SMDS_NodeIteratorPtr nodeIt = smDS->GetNodes();
1050 if ( nodeIt->more() )
1051 ngNodeId( nodeIt->next(), ngMesh, nodeNgIdMap );
1057 } // loop on submeshes
1060 nodeVec.resize( ngMesh.GetNP() + 1 );
1061 TNode2IdMap::iterator node_NgId, nodeNgIdEnd = nodeNgIdMap.end();
1062 for ( node_NgId = nodeNgIdMap.begin(); node_NgId != nodeNgIdEnd; ++node_NgId)
1063 nodeVec[ node_NgId->second ] = node_NgId->first;
1068 //================================================================================
1070 * \brief Duplicate mesh faces on internal geom faces
1072 //================================================================================
1074 void NETGENPlugin_Mesher::FixIntFaces(const netgen::OCCGeometry& occgeom,
1075 netgen::Mesh& ngMesh,
1076 NETGENPlugin_Internals& internalShapes)
1078 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1080 // find ng indices of internal faces
1082 for ( int ngFaceID = 1; ngFaceID <= occgeom.fmap.Extent(); ++ngFaceID )
1084 int smeshID = meshDS->ShapeToIndex( occgeom.fmap( ngFaceID ));
1085 if ( internalShapes.isInternalShape( smeshID ))
1086 ngFaceIds.insert( ngFaceID );
1088 if ( !ngFaceIds.empty() )
1091 int i, nbFaces = ngMesh.GetNSE();
1092 for (int i = 1; i <= nbFaces; ++i)
1094 netgen::Element2d elem = ngMesh.SurfaceElement(i);
1095 if ( ngFaceIds.count( elem.GetIndex() ))
1097 swap( elem[1], elem[2] );
1098 ngMesh.AddSurfaceElement (elem);
1106 //================================================================================
1107 // define gp_XY_Subtracted pointer to function calling gp_XY::Subtracted(gp_XY)
1108 gp_XY_FunPtr(Subtracted);
1109 //gp_XY_FunPtr(Added);
1111 //================================================================================
1113 * \brief Evaluate distance between two 2d points along the surface
1115 //================================================================================
1117 double evalDist( const gp_XY& uv1,
1119 const Handle(Geom_Surface)& surf,
1120 const int stopHandler=-1)
1122 if ( stopHandler > 0 ) // continue recursion
1124 gp_XY mid = SMESH_MesherHelper::GetMiddleUV( surf, uv1, uv2 );
1125 return evalDist( uv1,mid, surf, stopHandler-1 ) + evalDist( mid,uv2, surf, stopHandler-1 );
1127 double dist3D = surf->Value( uv1.X(), uv1.Y() ).Distance( surf->Value( uv2.X(), uv2.Y() ));
1128 if ( stopHandler == 0 ) // stop recursion
1131 // start recursion if necessary
1132 double dist2D = SMESH_MesherHelper::applyIn2D(surf, uv1, uv2, gp_XY_Subtracted, 0).Modulus();
1133 if ( fabs( dist3D - dist2D ) < dist2D * 1e-10 )
1134 return dist3D; // equal parametrization of a planar surface
1136 return evalDist( uv1, uv2, surf, 3 ); // start recursion
1139 //================================================================================
1141 * \brief Data of vertex internal in geom face
1143 //================================================================================
1147 gp_XY uv; //!< UV in face parametric space
1148 int ngId; //!< ng id of corrsponding node
1149 gp_XY uvClose; //!< UV of closest boundary node
1150 int ngIdClose; //!< ng id of closest boundary node
1153 //================================================================================
1155 * \brief Data of vertex internal in solid
1157 //================================================================================
1161 int ngId; //!< ng id of corresponding node
1162 int ngIdClose; //!< ng id of closest 2d mesh element
1163 int ngIdCloseN; //!< ng id of closest node of the closest 2d mesh element
1166 inline double dist2(const netgen::MeshPoint& p1, const netgen::MeshPoint& p2)
1168 return gp_Pnt( NGPOINT_COORDS(p1)).SquareDistance( gp_Pnt( NGPOINT_COORDS(p2)));
1172 //================================================================================
1174 * \brief Make netgen take internal vertices in faces into account by adding
1175 * segments including internal vertices
1177 * This function works in supposition that 1D mesh is already computed in ngMesh
1179 //================================================================================
1181 void NETGENPlugin_Mesher::AddIntVerticesInFaces(const netgen::OCCGeometry& occgeom,
1182 netgen::Mesh& ngMesh,
1183 vector<const SMDS_MeshNode*>& nodeVec,
1184 NETGENPlugin_Internals& internalShapes)
1186 if ( nodeVec.size() < ngMesh.GetNP() )
1187 nodeVec.resize( ngMesh.GetNP(), 0 );
1189 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1190 SMESH_MesherHelper helper( internalShapes.getMesh() );
1192 const map<int,list<int> >& face2Vert = internalShapes.getFacesWithVertices();
1193 map<int,list<int> >::const_iterator f2v = face2Vert.begin();
1194 for ( ; f2v != face2Vert.end(); ++f2v )
1196 const TopoDS_Face& face = TopoDS::Face( meshDS->IndexToShape( f2v->first ));
1197 if ( face.IsNull() ) continue;
1198 int faceNgID = occgeom.fmap.FindIndex (face);
1199 if ( faceNgID < 0 ) continue;
1201 TopLoc_Location loc;
1202 Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc);
1204 helper.SetSubShape( face );
1205 helper.SetElementsOnShape( true );
1207 // Get data of internal vertices and add them to ngMesh
1209 multimap< double, TIntVData > dist2VData; // sort vertices by distance from boundary nodes
1211 int i, nbSegInit = ngMesh.GetNSeg();
1213 // boundary characteristics
1214 double totSegLen2D = 0;
1217 const list<int>& iVertices = f2v->second;
1218 list<int>::const_iterator iv = iVertices.begin();
1219 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1222 // get node on vertex
1223 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1224 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1227 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1228 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1229 nV = SMESH_Algo::VertexNode( V, meshDS );
1230 if ( !nV ) continue;
1233 netgen::MeshPoint mp( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1234 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1235 vData.ngId = ngMesh.GetNP();
1236 nodeVec.push_back( nV );
1240 vData.uv = helper.GetNodeUV( face, nV, 0, &uvOK );
1241 if ( !uvOK ) helper.CheckNodeUV( face, nV, vData.uv, BRep_Tool::Tolerance(V),/*force=*/1);
1243 // loop on all segments of the face to find the node closest to vertex and to count
1244 // average segment 2d length
1245 double closeDist2 = numeric_limits<double>::max(), dist2;
1247 for (i = 1; i <= ngMesh.GetNSeg(); ++i)
1249 netgen::Segment & seg = ngMesh.LineSegment(i);
1250 if ( seg.si != faceNgID ) continue;
1252 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1254 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1255 if ( ngIdLast == seg[ iEnd ] ) continue;
1256 dist2 = helper.applyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1257 if ( dist2 < closeDist2 )
1258 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1259 ngIdLast = seg[ iEnd ];
1263 totSegLen2D += helper.applyIn2D(surf, uv[0], uv[1], gp_XY_Subtracted, false).Modulus();
1267 dist2VData.insert( make_pair( closeDist2, vData ));
1270 if ( totNbSeg == 0 ) break;
1271 double avgSegLen2d = totSegLen2D / totNbSeg;
1273 // Loop on vertices to add segments
1275 multimap< double, TIntVData >::iterator dist_vData = dist2VData.begin();
1276 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1278 double closeDist2 = dist_vData->first, dist2;
1279 TIntVData & vData = dist_vData->second;
1281 // try to find more close node among segments added for internal vertices
1282 for (i = nbSegInit+1; i <= ngMesh.GetNSeg(); ++i)
1284 netgen::Segment & seg = ngMesh.LineSegment(i);
1285 if ( seg.si != faceNgID ) continue;
1287 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1289 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1290 dist2 = helper.applyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1291 if ( dist2 < closeDist2 )
1292 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1295 // decide whether to use the closest node as the second end of segment or to
1296 // create a new point
1297 int segEnd1 = vData.ngId;
1298 int segEnd2 = vData.ngIdClose; // to use closest node
1299 gp_XY uvV = vData.uv, uvP = vData.uvClose;
1300 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1301 double nodeDist2D = sqrt( closeDist2 );
1302 double nodeDist3D = evalDist( vData.uv, vData.uvClose, surf );
1303 bool avgLenOK = ( avgSegLen2d < 0.75 * nodeDist2D );
1304 bool hintLenOK = ( segLenHint < 0.75 * nodeDist3D );
1305 //cout << "uvV " << uvV.X() <<","<<uvV.Y() << " ";
1306 if ( hintLenOK || avgLenOK )
1308 // create a point between the closest node and V
1311 double r = min( 0.5, ( hintLenOK ? segLenHint/nodeDist3D : avgSegLen2d/nodeDist2D ));
1312 // direction from V to closet node in 2D
1313 gp_Dir2d v2n( helper.applyIn2D(surf, uvP, uvV, gp_XY_Subtracted, false ));
1315 uvP = vData.uv + r * nodeDist2D * v2n.XY();
1316 gp_Pnt P = surf->Value( uvP.X(), uvP.Y() ).Transformed( loc );
1318 netgen::MeshPoint mp( netgen::Point<3> (P.X(), P.Y(), P.Z()));
1319 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1320 segEnd2 = ngMesh.GetNP();
1321 //cout << "Middle " << r << " uv " << uvP.X() << "," << uvP.Y() << "( " << ngMesh.Point(segEnd2).X()<<","<<ngMesh.Point(segEnd2).Y()<<","<<ngMesh.Point(segEnd2).Z()<<" )"<< endl;
1322 SMDS_MeshNode * nP = helper.AddNode(P.X(), P.Y(), P.Z());
1323 nodeVec.push_back( nP );
1325 //else cout << "at Node " << " uv " << uvP.X() << "," << uvP.Y() << endl;
1328 netgen::Segment seg;
1330 if ( segEnd1 > segEnd2 ) swap( segEnd1, segEnd2 ), swap( uvV, uvP );
1331 seg[0] = segEnd1; // ng node id
1332 seg[1] = segEnd2; // ng node id
1333 seg.edgenr = ngMesh.GetNSeg() + 1;// segment id
1336 seg.epgeominfo[ 0 ].dist = 0; // param on curve
1337 seg.epgeominfo[ 0 ].u = uvV.X();
1338 seg.epgeominfo[ 0 ].v = uvV.Y();
1339 seg.epgeominfo[ 1 ].dist = 1; // param on curve
1340 seg.epgeominfo[ 1 ].u = uvP.X();
1341 seg.epgeominfo[ 1 ].v = uvP.Y();
1343 // seg.epgeominfo[ 0 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1344 // seg.epgeominfo[ 1 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1346 ngMesh.AddSegment (seg);
1348 // add reverse segment
1349 swap (seg[0], seg[1]);
1350 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1351 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1352 ngMesh.AddSegment (seg);
1358 //================================================================================
1360 * \brief Make netgen take internal vertices in solids into account by adding
1361 * faces including internal vertices
1363 * This function works in supposition that 2D mesh is already computed in ngMesh
1365 //================================================================================
1367 void NETGENPlugin_Mesher::AddIntVerticesInSolids(const netgen::OCCGeometry& occgeom,
1368 netgen::Mesh& ngMesh,
1369 vector<const SMDS_MeshNode*>& nodeVec,
1370 NETGENPlugin_Internals& internalShapes)
1372 #ifdef DUMP_TRIANGLES_SCRIPT
1373 // create a python script making a mesh containing triangles added for internal vertices
1374 ofstream py(DUMP_TRIANGLES_SCRIPT);
1375 py << "import SMESH"<< endl
1376 << "from salome.smesh import smeshBuilder"<<endl
1377 << "smesh = smeshBuilder.New(salome.myStudy)"<<endl
1378 << "m = smesh.Mesh(name='triangles')" << endl;
1380 if ( nodeVec.size() < ngMesh.GetNP() )
1381 nodeVec.resize( ngMesh.GetNP(), 0 );
1383 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1384 SMESH_MesherHelper helper( internalShapes.getMesh() );
1386 const map<int,list<int> >& so2Vert = internalShapes.getSolidsWithVertices();
1387 map<int,list<int> >::const_iterator s2v = so2Vert.begin();
1388 for ( ; s2v != so2Vert.end(); ++s2v )
1390 const TopoDS_Shape& solid = meshDS->IndexToShape( s2v->first );
1391 if ( solid.IsNull() ) continue;
1392 int solidNgID = occgeom.somap.FindIndex (solid);
1393 if ( solidNgID < 0 && !occgeom.somap.IsEmpty() ) continue;
1395 helper.SetSubShape( solid );
1396 helper.SetElementsOnShape( true );
1398 // find ng indices of faces within the solid
1400 for (TopExp_Explorer fExp(solid, TopAbs_FACE); fExp.More(); fExp.Next() )
1401 ngFaceIds.insert( occgeom.fmap.FindIndex( fExp.Current() ));
1402 if ( ngFaceIds.size() == 1 && *ngFaceIds.begin() == 0 )
1403 ngFaceIds.insert( 1 );
1405 // Get data of internal vertices and add them to ngMesh
1407 multimap< double, TIntVSoData > dist2VData; // sort vertices by distance from ng faces
1409 int i, nbFaceInit = ngMesh.GetNSE();
1411 // boundary characteristics
1412 double totSegLen = 0;
1415 const list<int>& iVertices = s2v->second;
1416 list<int>::const_iterator iv = iVertices.begin();
1417 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1420 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1422 // get node on vertex
1423 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1426 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1427 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1428 nV = SMESH_Algo::VertexNode( V, meshDS );
1429 if ( !nV ) continue;
1432 netgen::MeshPoint mpV( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1433 ngMesh.AddPoint ( mpV, 1, netgen::FIXEDPOINT );
1434 vData.ngId = ngMesh.GetNP();
1435 nodeVec.push_back( nV );
1437 // loop on all 2d elements to find the one closest to vertex and to count
1438 // average segment length
1439 double closeDist2 = numeric_limits<double>::max(), avgDist2;
1440 for (i = 1; i <= ngMesh.GetNSE(); ++i)
1442 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1443 if ( !ngFaceIds.count( elem.GetIndex() )) continue;
1445 multimap< double, int> dist2nID; // sort nodes of element by distance from V
1446 for ( int j = 0; j < elem.GetNP(); ++j)
1448 netgen::MeshPoint mp = ngMesh.Point( elem[j] );
1449 double d2 = dist2( mpV, mp );
1450 dist2nID.insert( make_pair( d2, elem[j] ));
1451 avgDist2 += d2 / elem.GetNP();
1453 totNbSeg++, totSegLen+= sqrt( dist2( mp, ngMesh.Point( elem[(j+1)%elem.GetNP()])));
1455 double dist = dist2nID.begin()->first; //avgDist2;
1456 if ( dist < closeDist2 )
1457 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= dist;
1459 dist2VData.insert( make_pair( closeDist2, vData ));
1462 if ( totNbSeg == 0 ) break;
1463 double avgSegLen = totSegLen / totNbSeg;
1465 // Loop on vertices to add triangles
1467 multimap< double, TIntVSoData >::iterator dist_vData = dist2VData.begin();
1468 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1470 double closeDist2 = dist_vData->first;
1471 TIntVSoData & vData = dist_vData->second;
1473 const netgen::MeshPoint& mpV = ngMesh.Point( vData.ngId );
1475 // try to find more close face among ones added for internal vertices
1476 for (i = nbFaceInit+1; i <= ngMesh.GetNSE(); ++i)
1478 double avgDist2 = 0;
1479 multimap< double, int> dist2nID;
1480 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1481 for ( int j = 0; j < elem.GetNP(); ++j)
1483 double d = dist2( mpV, ngMesh.Point( elem[j] ));
1484 dist2nID.insert( make_pair( d, elem[j] ));
1485 avgDist2 += d / elem.GetNP();
1486 if ( avgDist2 < closeDist2 )
1487 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= avgDist2;
1490 // sort nodes of the closest face by angle with vector from V to the closest node
1491 const double tol = numeric_limits<double>::min();
1492 map< double, int > angle2ID;
1493 const netgen::Element2d& closeFace = ngMesh.SurfaceElement( vData.ngIdClose );
1494 netgen::MeshPoint mp[2];
1495 mp[0] = ngMesh.Point( vData.ngIdCloseN );
1496 gp_XYZ p1( NGPOINT_COORDS( mp[0] ));
1497 gp_XYZ pV( NGPOINT_COORDS( mpV ));
1498 gp_Vec v2p1( pV, p1 );
1499 double distN1 = v2p1.Magnitude();
1500 if ( distN1 <= tol ) continue;
1502 for ( int j = 0; j < closeFace.GetNP(); ++j)
1504 mp[1] = ngMesh.Point( closeFace[j] );
1505 gp_Vec v2p( pV, gp_Pnt( NGPOINT_COORDS( mp[1] )) );
1506 angle2ID.insert( make_pair( v2p1.Angle( v2p ), closeFace[j]));
1508 // get node with angle of 60 degrees or greater
1509 map< double, int >::iterator angle_id = angle2ID.lower_bound( 60. * M_PI / 180. );
1510 if ( angle_id == angle2ID.end() ) angle_id = --angle2ID.end();
1511 const double minAngle = 30. * M_PI / 180.;
1512 const double angle = angle_id->first;
1513 bool angleOK = ( angle > minAngle );
1515 // find points to create a triangle
1516 netgen::Element2d tri(3);
1518 tri[0] = vData.ngId;
1519 tri[1] = vData.ngIdCloseN; // to use the closest nodes
1520 tri[2] = angle_id->second; // to use the node with best angle
1522 // decide whether to use the closest node and the node with best angle or to create new ones
1523 for ( int isBestAngleN = 0; isBestAngleN < 2; ++isBestAngleN )
1525 bool createNew = !angleOK, distOK = true;
1527 int triInd = isBestAngleN ? 2 : 1;
1528 mp[isBestAngleN] = ngMesh.Point( tri[triInd] );
1533 double distN2 = sqrt( dist2( mpV, mp[isBestAngleN]));
1534 createNew = ( fabs( distN2 - distN1 ) > 0.25 * distN1 );
1536 else if ( angle < tol )
1538 v2p1.SetX( v2p1.X() + 1e-3 );
1544 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1545 bool avgLenOK = ( avgSegLen < 0.75 * distN1 );
1546 bool hintLenOK = ( segLenHint < 0.75 * distN1 );
1547 createNew = (createNew || avgLenOK || hintLenOK );
1548 // we create a new node not closer than 0.5 to the closest face
1549 // in order not to clash with other close face
1550 double r = min( 0.5, ( hintLenOK ? segLenHint : avgSegLen ) / distN1 );
1551 distFromV = r * distN1;
1555 // create a new point, between the node and the vertex if angleOK
1556 gp_XYZ p( NGPOINT_COORDS( mp[isBestAngleN] ));
1557 gp_Vec v2p( pV, p ); v2p.Normalize();
1558 if ( isBestAngleN && !angleOK )
1559 p = p1 + gp_Dir( v2p.XYZ() - v2p1.XYZ()).XYZ() * distN1 * 0.95;
1561 p = pV + v2p.XYZ() * distFromV;
1563 if ( !isBestAngleN ) p1 = p, distN1 = distFromV;
1565 mp[isBestAngleN].SetPoint( netgen::Point<3> (p.X(), p.Y(), p.Z()));
1566 ngMesh.AddPoint ( mp[isBestAngleN], 1, netgen::SURFACEPOINT );
1567 tri[triInd] = ngMesh.GetNP();
1568 nodeVec.push_back( helper.AddNode( p.X(), p.Y(), p.Z()) );
1571 ngMesh.AddSurfaceElement (tri);
1572 swap( tri[1], tri[2] );
1573 ngMesh.AddSurfaceElement (tri);
1575 #ifdef DUMP_TRIANGLES_SCRIPT
1576 py << "n1 = m.AddNode( "<< mpV(0)<<", "<< mpV(1)<<", "<< mpV(2)<<") "<< endl
1577 << "n2 = m.AddNode( "<< mp[0](0)<<", "<< mp[0](1)<<", "<< mp[0](2)<<") "<< endl
1578 << "n3 = m.AddNode( "<< mp[1](0)<<", "<< mp[1](1)<<", "<< mp[1](2)<<" )" << endl
1579 << "m.AddFace([n1,n2,n3])" << endl;
1581 } // loop on internal vertices of a solid
1583 } // loop on solids with internal vertices
1586 //================================================================================
1588 * \brief Fill netgen mesh with segments of a FACE
1589 * \param ngMesh - netgen mesh
1590 * \param geom - container of OCCT geometry to mesh
1591 * \param wires - data of nodes on FACE boundary
1592 * \param helper - mesher helper holding the FACE
1593 * \param nodeVec - vector of nodes in which node index == netgen ID
1594 * \retval SMESH_ComputeErrorPtr - error description
1596 //================================================================================
1598 SMESH_ComputeErrorPtr
1599 NETGENPlugin_Mesher::AddSegmentsToMesh(netgen::Mesh& ngMesh,
1600 netgen::OCCGeometry& geom,
1601 const TSideVector& wires,
1602 SMESH_MesherHelper& helper,
1603 vector< const SMDS_MeshNode* > & nodeVec,
1604 const bool overrideMinH)
1606 // ----------------------------
1607 // Check wires and count nodes
1608 // ----------------------------
1610 for ( int iW = 0; iW < wires.size(); ++iW )
1612 StdMeshers_FaceSidePtr wire = wires[ iW ];
1613 if ( wire->MissVertexNode() )
1615 // Commented for issue 0020960. It worked for the case, let's wait for case where it doesn't.
1616 // It seems that there is no reason for this limitation
1618 // (new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH, "Missing nodes on vertices"));
1620 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1621 if ( uvPtVec.size() != wire->NbPoints() )
1622 return SMESH_ComputeError::New(COMPERR_BAD_INPUT_MESH,
1623 SMESH_Comment("Unexpected nb of points on wire ") << iW
1624 << ": " << uvPtVec.size()<<" != "<<wire->NbPoints());
1625 nbNodes += wire->NbPoints();
1627 nodeVec.reserve( nodeVec.size() + nbNodes + 1 );
1628 if ( nodeVec.empty() )
1629 nodeVec.push_back( 0 );
1631 // -----------------
1633 // -----------------
1635 const bool wasNgMeshEmpty = ( ngMesh.GetNP() < 1 ); /* true => this method is called by
1636 NETGENPlugin_NETGEN_2D_ONLY */
1638 // map for nodes on vertices since they can be shared between wires
1639 // ( issue 0020676, face_int_box.brep) and nodes built by NETGEN
1640 map<const SMDS_MeshNode*, int > node2ngID;
1641 if ( !wasNgMeshEmpty ) // fill node2ngID with nodes built by NETGEN
1643 set< int > subIDs; // ids of sub-shapes of the FACE
1644 for ( int iW = 0; iW < wires.size(); ++iW )
1646 StdMeshers_FaceSidePtr wire = wires[ iW ];
1647 for ( int iE = 0, nbE = wire->NbEdges(); iE < nbE; ++iE )
1649 subIDs.insert( wire->EdgeID( iE ));
1650 subIDs.insert( helper.GetMeshDS()->ShapeToIndex( wire->FirstVertex( iE )));
1653 for ( size_t ngID = 1; ngID < nodeVec.size(); ++ngID )
1654 if ( subIDs.count( nodeVec[ngID]->getshapeId() ))
1655 node2ngID.insert( make_pair( nodeVec[ngID], ngID ));
1658 const int solidID = 0, faceID = geom.fmap.FindIndex( helper.GetSubShape() );
1659 if ( ngMesh.GetNFD() < 1 )
1660 ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(faceID, solidID, solidID, 0));
1662 for ( int iW = 0; iW < wires.size(); ++iW )
1664 StdMeshers_FaceSidePtr wire = wires[ iW ];
1665 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1666 const int nbSegments = wire->NbPoints() - 1;
1668 // assure the 1st node to be in node2ngID, which is needed to correctly
1669 // "close chain of segments" (see below) in case if the 1st node is not
1670 // onVertex because it is on a Viscous layer
1671 node2ngID.insert( make_pair( uvPtVec[ 0 ].node, ngMesh.GetNP() + 1 ));
1673 // compute length of every segment
1674 vector<double> segLen( nbSegments );
1675 for ( int i = 0; i < nbSegments; ++i )
1676 segLen[i] = SMESH_TNodeXYZ( uvPtVec[ i ].node ).Distance( uvPtVec[ i+1 ].node );
1678 int edgeID = 1, posID = -2;
1679 bool isInternalWire = false;
1680 double vertexNormPar = 0;
1681 const int prevNbNGSeg = ngMesh.GetNSeg();
1682 for ( int i = 0; i < nbSegments; ++i ) // loop on segments
1684 // Add the first point of a segment
1686 const SMDS_MeshNode * n = uvPtVec[ i ].node;
1687 const int posShapeID = n->getshapeId();
1688 bool onVertex = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX );
1689 bool onEdge = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_EDGE );
1691 // skip nodes on degenerated edges
1692 if ( helper.IsDegenShape( posShapeID ) &&
1693 helper.IsDegenShape( uvPtVec[ i+1 ].node->getshapeId() ))
1696 int ngID1 = ngMesh.GetNP() + 1, ngID2 = ngID1+1;
1697 if ( onVertex || ( !wasNgMeshEmpty && onEdge ))
1698 ngID1 = node2ngID.insert( make_pair( n, ngID1 )).first->second;
1699 if ( ngID1 > ngMesh.GetNP() )
1701 netgen::MeshPoint mp( netgen::Point<3> (n->X(), n->Y(), n->Z()) );
1702 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1703 nodeVec.push_back( n );
1705 else // n is in ngMesh already, and ngID2 in prev segment is wrong
1707 ngID2 = ngMesh.GetNP() + 1;
1708 if ( i > 0 ) // prev segment belongs to same wire
1710 netgen::Segment& prevSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
1717 netgen::Segment seg;
1719 seg[0] = ngID1; // ng node id
1720 seg[1] = ngID2; // ng node id
1721 seg.edgenr = ngMesh.GetNSeg() + 1; // ng segment id
1722 seg.si = faceID; // = geom.fmap.FindIndex (face);
1724 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1726 const UVPtStruct& pnt = uvPtVec[ i + iEnd ];
1728 seg.epgeominfo[ iEnd ].dist = pnt.param; // param on curve
1729 seg.epgeominfo[ iEnd ].u = pnt.u;
1730 seg.epgeominfo[ iEnd ].v = pnt.v;
1732 // find out edge id and node parameter on edge
1733 onVertex = ( pnt.normParam + 1e-10 > vertexNormPar );
1734 if ( onVertex || posShapeID != posID )
1737 double normParam = pnt.normParam;
1739 normParam = 0.5 * ( uvPtVec[ i ].normParam + uvPtVec[ i+1 ].normParam );
1740 int edgeIndexInWire = wire->EdgeIndex( normParam );
1741 vertexNormPar = wire->LastParameter( edgeIndexInWire );
1742 const TopoDS_Edge& edge = wire->Edge( edgeIndexInWire );
1743 edgeID = geom.emap.FindIndex( edge );
1745 isInternalWire = ( edge.Orientation() == TopAbs_INTERNAL );
1746 // if ( onVertex ) // param on curve is different on each of two edges
1747 // seg.epgeominfo[ iEnd ].dist = helper.GetNodeU( edge, pnt.node );
1749 seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
1752 ngMesh.AddSegment (seg);
1754 // restrict size of elements near the segment
1755 SMESH_TNodeXYZ np1( n ), np2( uvPtVec[ i+1 ].node );
1756 // get an average size of adjacent segments to avoid sharp change of
1757 // element size (regression on issue 0020452, note 0010898)
1758 int iPrev = SMESH_MesherHelper::WrapIndex( i-1, nbSegments );
1759 int iNext = SMESH_MesherHelper::WrapIndex( i+1, nbSegments );
1760 double sumH = segLen[ iPrev ] + segLen[ i ] + segLen[ iNext ];
1761 int nbSeg = ( int( segLen[ iPrev ] > sumH / 100.) +
1762 int( segLen[ i ] > sumH / 100.) +
1763 int( segLen[ iNext ] > sumH / 100.));
1765 RestrictLocalSize( ngMesh, 0.5*(np1+np2), sumH / nbSeg, overrideMinH );
1767 if ( isInternalWire )
1769 swap (seg[0], seg[1]);
1770 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1771 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1772 ngMesh.AddSegment (seg);
1774 } // loop on segments on a wire
1776 // close chain of segments
1777 if ( nbSegments > 0 )
1779 netgen::Segment& lastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() - int( isInternalWire));
1780 const SMDS_MeshNode * lastNode = uvPtVec.back().node;
1781 lastSeg[1] = node2ngID.insert( make_pair( lastNode, lastSeg[1] )).first->second;
1782 if ( lastSeg[1] > ngMesh.GetNP() )
1784 netgen::MeshPoint mp( netgen::Point<3> (lastNode->X(), lastNode->Y(), lastNode->Z()) );
1785 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1786 nodeVec.push_back( lastNode );
1788 if ( isInternalWire )
1790 netgen::Segment& realLastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
1791 realLastSeg[0] = lastSeg[1];
1795 #ifdef DUMP_SEGMENTS
1796 cout << "BEGIN WIRE " << iW << endl;
1797 for ( int i = prevNbNGSeg+1; i <= ngMesh.GetNSeg(); ++i )
1799 netgen::Segment& seg = ngMesh.LineSegment( i );
1801 netgen::Segment& prevSeg = ngMesh.LineSegment( i-1 );
1802 if ( seg[0] == prevSeg[1] && seg[1] == prevSeg[0] )
1804 cout << "Segment: " << seg.edgenr << endl << "\tis REVRESE of the previous one" << endl;
1808 cout << "Segment: " << seg.edgenr << endl
1809 << "\tp1: " << seg[0] << endl
1810 << "\tp2: " << seg[1] << endl
1811 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
1812 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
1813 << "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
1814 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
1815 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl
1816 << "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
1818 cout << "--END WIRE " << iW << endl;
1821 } // loop on WIREs of a FACE
1823 // add a segment instead of an internal vertex
1824 if ( wasNgMeshEmpty )
1826 NETGENPlugin_Internals intShapes( *helper.GetMesh(), helper.GetSubShape(), /*is3D=*/false );
1827 AddIntVerticesInFaces( geom, ngMesh, nodeVec, intShapes );
1829 ngMesh.CalcSurfacesOfNode();
1834 //================================================================================
1836 * \brief Fill SMESH mesh according to contents of netgen mesh
1837 * \param occgeo - container of OCCT geometry to mesh
1838 * \param ngMesh - netgen mesh
1839 * \param initState - bn of entities in netgen mesh before computing
1840 * \param sMesh - SMESH mesh to fill in
1841 * \param nodeVec - vector of nodes in which node index == netgen ID
1842 * \param comment - returns problem description
1843 * \param quadHelper - holder of medium nodes of sub-meshes
1844 * \retval int - error
1846 //================================================================================
1848 int NETGENPlugin_Mesher::FillSMesh(const netgen::OCCGeometry& occgeo,
1849 netgen::Mesh& ngMesh,
1850 const NETGENPlugin_ngMeshInfo& initState,
1852 std::vector<const SMDS_MeshNode*>& nodeVec,
1853 SMESH_Comment& comment,
1854 SMESH_MesherHelper* quadHelper)
1856 int nbNod = ngMesh.GetNP();
1857 int nbSeg = ngMesh.GetNSeg();
1858 int nbFac = ngMesh.GetNSE();
1859 int nbVol = ngMesh.GetNE();
1861 SMESHDS_Mesh* meshDS = sMesh.GetMeshDS();
1863 // quadHelper is used for either
1864 // 1) making quadratic elements when a lower dimention mesh is loaded
1865 // to SMESH before convertion to quadratic by NETGEN
1866 // 2) sewing of quadratic elements with quadratic elements of sub-meshes
1867 if ( quadHelper && !quadHelper->GetIsQuadratic() && quadHelper->GetTLinkNodeMap().empty() )
1870 // -------------------------------------
1871 // Create and insert nodes into nodeVec
1872 // -------------------------------------
1874 nodeVec.resize( nbNod + 1 );
1875 int i, nbInitNod = initState._nbNodes;
1876 for (i = nbInitNod+1; i <= nbNod; ++i )
1878 const netgen::MeshPoint& ngPoint = ngMesh.Point(i);
1879 SMDS_MeshNode* node = NULL;
1880 TopoDS_Vertex aVert;
1881 // First, netgen creates nodes on vertices in occgeo.vmap,
1882 // so node index corresponds to vertex index
1883 // but (issue 0020776) netgen does not create nodes with equal coordinates
1884 if ( i-nbInitNod <= occgeo.vmap.Extent() )
1886 gp_Pnt p ( NGPOINT_COORDS(ngPoint) );
1887 for (int iV = i-nbInitNod; aVert.IsNull() && iV <= occgeo.vmap.Extent(); ++iV)
1889 aVert = TopoDS::Vertex( occgeo.vmap( iV ));
1890 gp_Pnt pV = BRep_Tool::Pnt( aVert );
1891 if ( p.SquareDistance( pV ) > 1e-20 )
1894 node = const_cast<SMDS_MeshNode*>( SMESH_Algo::VertexNode( aVert, meshDS ));
1897 if (!node) // node not found on vertex
1899 node = meshDS->AddNode( NGPOINT_COORDS( ngPoint ));
1900 if (!aVert.IsNull())
1901 meshDS->SetNodeOnVertex(node, aVert);
1906 // -------------------------------------------
1907 // Create mesh segments along geometric edges
1908 // -------------------------------------------
1910 int nbInitSeg = initState._nbSegments;
1911 for (i = nbInitSeg+1; i <= nbSeg; ++i )
1913 const netgen::Segment& seg = ngMesh.LineSegment(i);
1915 int pinds[3] = { seg.pnums[0], seg.pnums[1], seg.pnums[2] };
1918 for (int j=0; j < 3; ++j)
1920 int pind = pinds[j];
1921 if (pind <= 0 || !nodeVec_ACCESS(pind))
1929 int aGeomEdgeInd = seg.epgeominfo[j].edgenr;
1930 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
1931 aEdge = TopoDS::Edge(occgeo.emap(aGeomEdgeInd));
1933 param = seg.epgeominfo[j].dist;
1936 else // middle point
1938 param = param2 * 0.5;
1940 if (!aEdge.IsNull() && nodeVec_ACCESS(pind)->getshapeId() < 1)
1942 meshDS->SetNodeOnEdge(nodeVec_ACCESS(pind), aEdge, param);
1947 SMDS_MeshEdge* edge = 0;
1948 if (nbp == 2) // second order ?
1950 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1])))
1952 if ( quadHelper ) // final mesh must be quadratic
1953 edge = quadHelper->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
1955 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
1959 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
1960 nodeVec_ACCESS(pinds[2])))
1962 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
1963 nodeVec_ACCESS(pinds[2]));
1967 if ( comment.empty() ) comment << "Cannot create a mesh edge";
1968 MESSAGE("Cannot create a mesh edge");
1969 nbSeg = nbFac = nbVol = 0;
1972 if ( !aEdge.IsNull() && edge->getshapeId() < 1 )
1973 meshDS->SetMeshElementOnShape(edge, aEdge);
1975 else if ( comment.empty() )
1977 comment << "Invalid netgen segment #" << i;
1981 // ----------------------------------------
1982 // Create mesh faces along geometric faces
1983 // ----------------------------------------
1985 int nbInitFac = initState._nbFaces;
1986 int quadFaceID = ngMesh.GetNFD() + 1;
1987 if ( nbInitFac < nbFac )
1988 // add a faces descriptor to exclude qudrangle elements generated by NETGEN
1989 // from computation of 3D mesh
1990 ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(quadFaceID, /*solid1=*/0, /*solid2=*/0, 0));
1992 vector<const SMDS_MeshNode*> nodes;
1993 for (i = nbInitFac+1; i <= nbFac; ++i )
1995 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1996 int aGeomFaceInd = elem.GetIndex();
1998 if (aGeomFaceInd > 0 && aGeomFaceInd <= occgeo.fmap.Extent())
1999 aFace = TopoDS::Face(occgeo.fmap(aGeomFaceInd));
2001 for (int j=1; j <= elem.GetNP(); ++j)
2003 int pind = elem.PNum(j);
2004 if ( pind < 1 || pind >= nodeVec.size() )
2006 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind))
2008 nodes.push_back( node );
2009 if (!aFace.IsNull() && node->getshapeId() < 1)
2011 const netgen::PointGeomInfo& pgi = elem.GeomInfoPi(j);
2012 meshDS->SetNodeOnFace(node, aFace, pgi.u, pgi.v);
2016 if ( nodes.size() != elem.GetNP() )
2018 if ( comment.empty() )
2019 comment << "Invalid netgen 2d element #" << i;
2020 continue; // bad node ids
2022 SMDS_MeshFace* face = NULL;
2023 switch (elem.GetType())
2026 if ( quadHelper ) // final mesh must be quadratic
2027 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2]);
2029 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2]);
2032 if ( quadHelper ) // final mesh must be quadratic
2033 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2035 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2036 // exclude qudrangle elements from computation of 3D mesh
2037 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2040 nodes[5] = mediumNode( nodes[0],nodes[1],nodes[5], quadHelper );
2041 nodes[3] = mediumNode( nodes[1],nodes[2],nodes[3], quadHelper );
2042 nodes[4] = mediumNode( nodes[2],nodes[0],nodes[4], quadHelper );
2043 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[5],nodes[3],nodes[4]);
2046 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2047 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2048 nodes[5] = mediumNode( nodes[2],nodes[3],nodes[5], quadHelper );
2049 nodes[6] = mediumNode( nodes[3],nodes[0],nodes[6], quadHelper );
2050 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3],
2051 nodes[4],nodes[7],nodes[5],nodes[6]);
2052 // exclude qudrangle elements from computation of 3D mesh
2053 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2056 MESSAGE("NETGEN created a face of unexpected type, ignoring");
2061 if ( comment.empty() ) comment << "Cannot create a mesh face";
2062 MESSAGE("Cannot create a mesh face");
2063 nbSeg = nbFac = nbVol = 0;
2066 if (!aFace.IsNull())
2067 meshDS->SetMeshElementOnShape(face, aFace);
2070 // ------------------
2071 // Create tetrahedra
2072 // ------------------
2074 for (i = 1; i <= nbVol; ++i)
2076 const netgen::Element& elem = ngMesh.VolumeElement(i);
2077 int aSolidInd = elem.GetIndex();
2078 TopoDS_Solid aSolid;
2079 if (aSolidInd > 0 && aSolidInd <= occgeo.somap.Extent())
2080 aSolid = TopoDS::Solid(occgeo.somap(aSolidInd));
2082 for (int j=1; j <= elem.GetNP(); ++j)
2084 int pind = elem.PNum(j);
2085 if ( pind < 1 || pind >= nodeVec.size() )
2087 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind) )
2089 nodes.push_back(node);
2090 if ( !aSolid.IsNull() && node->getshapeId() < 1 )
2091 meshDS->SetNodeInVolume(node, aSolid);
2094 if ( nodes.size() != elem.GetNP() )
2096 if ( comment.empty() )
2097 comment << "Invalid netgen 3d element #" << i;
2100 SMDS_MeshVolume* vol = NULL;
2101 switch (elem.GetType())
2104 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3]);
2107 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2108 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2109 nodes[5] = mediumNode( nodes[2],nodes[0],nodes[5], quadHelper );
2110 nodes[6] = mediumNode( nodes[0],nodes[3],nodes[6], quadHelper );
2111 nodes[8] = mediumNode( nodes[1],nodes[3],nodes[8], quadHelper );
2112 nodes[9] = mediumNode( nodes[2],nodes[3],nodes[9], quadHelper );
2113 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3],
2114 nodes[4],nodes[7],nodes[5],nodes[6],nodes[8],nodes[9]);
2117 MESSAGE("NETGEN created a volume of unexpected type, ignoring");
2122 if ( comment.empty() ) comment << "Cannot create a mesh volume";
2123 MESSAGE("Cannot create a mesh volume");
2124 nbSeg = nbFac = nbVol = 0;
2127 if (!aSolid.IsNull())
2128 meshDS->SetMeshElementOnShape(vol, aSolid);
2130 return comment.empty() ? 0 : 1;
2135 //================================================================================
2137 * \brief Restrict size of elements on the given edge
2139 //================================================================================
2141 void setLocalSize(const TopoDS_Edge& edge,
2145 const int nb = 1000;
2146 Standard_Real u1, u2;
2147 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, u1, u2);
2148 if ( curve.IsNull() )
2150 TopoDS_Iterator vIt( edge );
2151 if ( !vIt.More() ) return;
2152 gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex( vIt.Value() ));
2153 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
2157 Standard_Real delta = (u2-u1)/nb;
2158 for(int i=0; i<nb; i++)
2160 Standard_Real u = u1 + delta*i;
2161 gp_Pnt p = curve->Value(u);
2162 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
2163 netgen::Point3d pi(p.X(), p.Y(), p.Z());
2164 double resultSize = mesh.GetH(pi);
2165 if ( resultSize - size > 0.1*size )
2166 // netgen does restriction iff oldH/newH > 1.2 (localh.cpp:136)
2167 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), resultSize/1.201 );
2172 //================================================================================
2174 * \brief Convert error into text
2176 //================================================================================
2178 std::string text(int err)
2183 SMESH_Comment("Error in netgen::OCCGenerateMesh() at ") << netgen::multithread.task;
2186 //================================================================================
2188 * \brief Convert exception into text
2190 //================================================================================
2192 std::string text(Standard_Failure& ex)
2194 SMESH_Comment str("Exception in netgen::OCCGenerateMesh()");
2195 str << " at " << netgen::multithread.task
2196 << ": " << ex.DynamicType()->Name();
2197 if ( ex.GetMessageString() && strlen( ex.GetMessageString() ))
2198 str << ": " << ex.GetMessageString();
2201 //================================================================================
2203 * \brief Convert exception into text
2205 //================================================================================
2207 std::string text(netgen::NgException& ex)
2209 SMESH_Comment str("NgException");
2210 if ( strlen( netgen::multithread.task ) > 0 )
2211 str << " at " << netgen::multithread.task;
2212 str << ": " << ex.What();
2216 //================================================================================
2218 * \brief Looks for triangles lying on a SOLID
2220 //================================================================================
2222 bool hasBadElemOnSolid( const list<const SMDS_MeshElement*>& elems,
2223 SMESH_subMesh* solidSM )
2225 TopTools_IndexedMapOfShape solidSubs;
2226 TopExp::MapShapes( solidSM->GetSubShape(), solidSubs );
2227 SMESHDS_Mesh* mesh = solidSM->GetFather()->GetMeshDS();
2229 list<const SMDS_MeshElement*>::const_iterator e = elems.begin();
2230 for ( ; e != elems.end(); ++e )
2232 const SMDS_MeshElement* elem = *e;
2233 // if ( elem->GetType() != SMDSAbs_Face ) -- 23047
2235 int nbNodesOnSolid = 0, nbNodes = elem->NbNodes();
2236 SMDS_NodeIteratorPtr nIt = elem->nodeIterator();
2237 while ( nIt->more() )
2239 const SMDS_MeshNode* n = nIt->next();
2240 const TopoDS_Shape& s = mesh->IndexToShape( n->getshapeId() );
2241 nbNodesOnSolid += ( !s.IsNull() && solidSubs.Contains( s ));
2242 if ( nbNodesOnSolid > 2 ||
2243 nbNodesOnSolid == nbNodes)
2250 const double edgeMeshingTime = 0.001;
2251 const double faceMeshingTime = 0.019;
2252 const double edgeFaceMeshingTime = edgeMeshingTime + faceMeshingTime;
2253 const double faceOptimizTime = 0.06;
2254 const double voluMeshingTime = 0.15;
2255 const double volOptimizeTime = 0.77;
2258 //=============================================================================
2260 * Here we are going to use the NETGEN mesher
2262 //=============================================================================
2264 bool NETGENPlugin_Mesher::Compute()
2266 NETGENPlugin_NetgenLibWrapper ngLib;
2268 netgen::MeshingParameters& mparams = netgen::mparam;
2269 MESSAGE("Compute with:\n"
2270 " max size = " << mparams.maxh << "\n"
2271 " segments per edge = " << mparams.segmentsperedge);
2273 " growth rate = " << mparams.grading << "\n"
2274 " elements per radius = " << mparams.curvaturesafety << "\n"
2275 " second order = " << mparams.secondorder << "\n"
2276 " quad allowed = " << mparams.quad << "\n"
2277 " surface curvature = " << mparams.uselocalh << "\n"
2278 " fuse edges = " << netgen::merge_solids);
2280 SMESH_ComputeErrorPtr error = SMESH_ComputeError::New();
2281 SMESH_MesherHelper quadHelper( *_mesh );
2282 quadHelper.SetIsQuadratic( mparams.secondorder );
2284 static string debugFile = "/tmp/ngMesh.py"; /* to call toPython( _ngMesh, debugFile )
2285 while debugging netgen */
2286 // -------------------------
2287 // Prepare OCC geometry
2288 // -------------------------
2290 netgen::OCCGeometry occgeo;
2291 list< SMESH_subMesh* > meshedSM[3]; // for 0-2 dimensions
2292 NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
2293 PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
2296 _totalTime = edgeFaceMeshingTime;
2298 _totalTime += faceOptimizTime;
2300 _totalTime += voluMeshingTime + ( _optimize ? volOptimizeTime : 0 );
2301 double doneTime = 0;
2304 _curShapeIndex = -1;
2306 // -------------------------
2307 // Generate the mesh
2308 // -------------------------
2311 NETGENPlugin_ngMeshInfo initState; // it remembers size of ng mesh equal to size of Smesh
2313 SMESH_Comment comment;
2316 // vector of nodes in which node index == netgen ID
2317 vector< const SMDS_MeshNode* > nodeVec;
2325 // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
2326 mparams.uselocalh = false;
2327 mparams.grading = 0.8; // not limitited size growth
2329 if ( _simpleHyp->GetNumberOfSegments() )
2331 mparams.maxh = occgeo.boundingbox.Diam();
2334 mparams.maxh = _simpleHyp->GetLocalLength();
2337 if ( mparams.maxh == 0.0 )
2338 mparams.maxh = occgeo.boundingbox.Diam();
2339 if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
2340 mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
2342 // Local size on faces
2343 occgeo.face_maxh = mparams.maxh;
2345 // Let netgen create _ngMesh and calculate element size on not meshed shapes
2349 int startWith = netgen::MESHCONST_ANALYSE;
2350 int endWith = netgen::MESHCONST_ANALYSE;
2355 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2357 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2359 if(netgen::multithread.terminate)
2362 comment << text(err);
2364 catch (Standard_Failure& ex)
2366 comment << text(ex);
2368 err = 0; //- MESHCONST_ANALYSE isn't so important step
2371 ngLib.setMesh(( Ng_Mesh*) _ngMesh );
2373 _ngMesh->ClearFaceDescriptors(); // we make descriptors our-self
2377 // Pass 1D simple parameters to NETGEN
2378 // --------------------------------
2379 int nbSeg = _simpleHyp->GetNumberOfSegments();
2380 double segSize = _simpleHyp->GetLocalLength();
2381 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
2383 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
2385 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
2386 setLocalSize( e, segSize, *_ngMesh );
2389 else // if ( ! _simpleHyp )
2391 // Local size on vertices and edges
2392 // --------------------------------
2393 for(std::map<int,double>::const_iterator it=EdgeId2LocalSize.begin(); it!=EdgeId2LocalSize.end(); it++)
2395 int key = (*it).first;
2396 double hi = (*it).second;
2397 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2398 const TopoDS_Edge& e = TopoDS::Edge(shape);
2399 setLocalSize( e, hi, *_ngMesh );
2401 for(std::map<int,double>::const_iterator it=VertexId2LocalSize.begin(); it!=VertexId2LocalSize.end(); it++)
2403 int key = (*it).first;
2404 double hi = (*it).second;
2405 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2406 const TopoDS_Vertex& v = TopoDS::Vertex(shape);
2407 gp_Pnt p = BRep_Tool::Pnt(v);
2408 NETGENPlugin_Mesher::RestrictLocalSize( *_ngMesh, p.XYZ(), hi );
2410 for(map<int,double>::const_iterator it=FaceId2LocalSize.begin();
2411 it!=FaceId2LocalSize.end(); it++)
2413 int key = (*it).first;
2414 double val = (*it).second;
2415 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2416 int faceNgID = occgeo.fmap.FindIndex(shape);
2417 occgeo.SetFaceMaxH(faceNgID, val);
2418 for ( TopExp_Explorer edgeExp( shape, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
2419 setLocalSize( TopoDS::Edge( edgeExp.Current() ), val, *_ngMesh );
2423 // Precompute internal edges (issue 0020676) in order to
2424 // add mesh on them correctly (twice) to netgen mesh
2425 if ( !err && internals.hasInternalEdges() )
2427 // load internal shapes into OCCGeometry
2428 netgen::OCCGeometry intOccgeo;
2429 internals.getInternalEdges( intOccgeo.fmap, intOccgeo.emap, intOccgeo.vmap, meshedSM );
2430 intOccgeo.boundingbox = occgeo.boundingbox;
2431 intOccgeo.shape = occgeo.shape;
2432 intOccgeo.face_maxh.SetSize(intOccgeo.fmap.Extent());
2433 intOccgeo.face_maxh = netgen::mparam.maxh;
2434 netgen::Mesh *tmpNgMesh = NULL;
2438 // compute local H on internal shapes in the main mesh
2439 //OCCSetLocalMeshSize(intOccgeo, *_ngMesh); it deletes _ngMesh->localH
2441 // let netgen create a temporary mesh
2443 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2445 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2447 if(netgen::multithread.terminate)
2450 // copy LocalH from the main to temporary mesh
2451 initState.transferLocalH( _ngMesh, tmpNgMesh );
2453 // compute mesh on internal edges
2454 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2456 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2458 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2460 comment << text(err);
2462 catch (Standard_Failure& ex)
2464 comment << text(ex);
2467 initState.restoreLocalH( tmpNgMesh );
2469 // fill SMESH by netgen mesh
2470 vector< const SMDS_MeshNode* > tmpNodeVec;
2471 FillSMesh( intOccgeo, *tmpNgMesh, initState, *_mesh, tmpNodeVec, comment );
2472 err = ( err || !comment.empty() );
2474 nglib::Ng_DeleteMesh((nglib::Ng_Mesh*)tmpNgMesh);
2477 // Fill _ngMesh with nodes and segments of computed submeshes
2480 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_0D ]) &&
2481 FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_1D ], &quadHelper));
2483 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2488 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2493 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2495 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2497 if(netgen::multithread.terminate)
2500 comment << text(err);
2502 catch (Standard_Failure& ex)
2504 comment << text(ex);
2509 _ticTime = ( doneTime += edgeMeshingTime ) / _totalTime / _progressTic;
2511 mparams.uselocalh = true; // restore as it is used at surface optimization
2513 // ---------------------
2514 // compute surface mesh
2515 // ---------------------
2518 // Pass 2D simple parameters to NETGEN
2520 if ( double area = _simpleHyp->GetMaxElementArea() ) {
2522 mparams.maxh = sqrt(2. * area/sqrt(3.0));
2523 mparams.grading = 0.4; // moderate size growth
2526 // length from edges
2527 if ( _ngMesh->GetNSeg() ) {
2528 double edgeLength = 0;
2529 TopTools_MapOfShape visitedEdges;
2530 for ( TopExp_Explorer exp( _shape, TopAbs_EDGE ); exp.More(); exp.Next() )
2531 if( visitedEdges.Add(exp.Current()) )
2532 edgeLength += SMESH_Algo::EdgeLength( TopoDS::Edge( exp.Current() ));
2533 // we have to multiply length by 2 since for each TopoDS_Edge there
2534 // are double set of NETGEN edges, in other words, we have to
2535 // divide _ngMesh->GetNSeg() by 2.
2536 mparams.maxh = 2*edgeLength / _ngMesh->GetNSeg();
2539 mparams.maxh = 1000;
2541 mparams.grading = 0.2; // slow size growth
2543 mparams.quad = _simpleHyp->GetAllowQuadrangles();
2544 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2545 _ngMesh->SetGlobalH (mparams.maxh);
2546 netgen::Box<3> bb = occgeo.GetBoundingBox();
2547 bb.Increase (bb.Diam()/20);
2548 _ngMesh->SetLocalH (bb.PMin(), bb.PMax(), mparams.grading);
2551 // Care of vertices internal in faces (issue 0020676)
2552 if ( internals.hasInternalVertexInFace() )
2554 // store computed segments in SMESH in order not to create SMESH
2555 // edges for ng segments added by AddIntVerticesInFaces()
2556 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2557 // add segments to faces with internal vertices
2558 AddIntVerticesInFaces( occgeo, *_ngMesh, nodeVec, internals );
2559 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2562 // Build viscous layers
2563 if ( _isViscousLayers2D )
2565 if ( !internals.hasInternalVertexInFace() ) {
2566 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2567 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2569 SMESH_ProxyMesh::Ptr viscousMesh;
2570 SMESH_MesherHelper helper( *_mesh );
2571 for ( int faceID = 1; faceID <= occgeo.fmap.Extent(); ++faceID )
2573 const TopoDS_Face& F = TopoDS::Face( occgeo.fmap( faceID ));
2574 viscousMesh = StdMeshers_ViscousLayers2D::Compute( *_mesh, F );
2577 // exclude from computation ng segments built on EDGEs of F
2578 for (int i = 1; i <= _ngMesh->GetNSeg(); i++)
2580 netgen::Segment & seg = _ngMesh->LineSegment(i);
2581 if (seg.si == faceID)
2584 // add new segments to _ngMesh instead of excluded ones
2585 helper.SetSubShape( F );
2587 StdMeshers_FaceSide::GetFaceWires( F, *_mesh, /*skipMediumNodes=*/true,
2588 error, viscousMesh );
2589 error = AddSegmentsToMesh( *_ngMesh, occgeo, wires, helper, nodeVec );
2591 if ( !error ) error = SMESH_ComputeError::New();
2593 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2596 // Let netgen compute 2D mesh
2597 startWith = netgen::MESHCONST_MESHSURFACE;
2598 endWith = _optimize ? netgen::MESHCONST_OPTSURFACE : netgen::MESHCONST_MESHSURFACE;
2603 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2605 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2607 if(netgen::multithread.terminate)
2610 comment << text (err);
2612 catch (Standard_Failure& ex)
2614 comment << text(ex);
2615 //err = 1; -- try to make volumes anyway
2617 catch (netgen::NgException exc)
2619 comment << text(exc);
2620 //err = 1; -- try to make volumes anyway
2625 doneTime += faceMeshingTime + ( _optimize ? faceOptimizTime : 0 );
2626 _ticTime = doneTime / _totalTime / _progressTic;
2628 // ---------------------
2629 // generate volume mesh
2630 // ---------------------
2631 // Fill _ngMesh with nodes and faces of computed 2D submeshes
2632 if ( !err && _isVolume && ( !meshedSM[ MeshDim_2D ].empty() || mparams.quad ))
2634 // load SMESH with computed segments and faces
2635 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2637 // compute pyramids on quadrangles
2638 SMESH_ProxyMesh::Ptr proxyMesh;
2639 if ( _mesh->NbQuadrangles() > 0 )
2640 for ( int iS = 1; iS <= occgeo.somap.Extent(); ++iS )
2642 StdMeshers_QuadToTriaAdaptor* Adaptor = new StdMeshers_QuadToTriaAdaptor;
2643 proxyMesh.reset( Adaptor );
2645 int nbPyrams = _mesh->NbPyramids();
2646 Adaptor->Compute( *_mesh, occgeo.somap(iS) );
2647 if ( nbPyrams != _mesh->NbPyramids() )
2649 list< SMESH_subMesh* > quadFaceSM;
2650 for (TopExp_Explorer face(occgeo.somap(iS), TopAbs_FACE); face.More(); face.Next())
2651 if ( Adaptor->GetProxySubMesh( face.Current() ))
2653 quadFaceSM.push_back( _mesh->GetSubMesh( face.Current() ));
2654 meshedSM[ MeshDim_2D ].remove( quadFaceSM.back() );
2656 FillNgMesh(occgeo, *_ngMesh, nodeVec, quadFaceSM, &quadHelper, proxyMesh);
2659 // fill _ngMesh with faces of sub-meshes
2660 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_2D ], &quadHelper));
2661 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2662 //toPython( _ngMesh, "/tmp/ngPython.py");
2664 if (!err && _isVolume)
2666 // Pass 3D simple parameters to NETGEN
2667 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
2668 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
2670 if ( double vol = simple3d->GetMaxElementVolume() ) {
2672 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
2673 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2676 // length from faces
2677 mparams.maxh = _ngMesh->AverageH();
2679 _ngMesh->SetGlobalH (mparams.maxh);
2680 mparams.grading = 0.4;
2682 _ngMesh->CalcLocalH(mparams.grading);
2684 _ngMesh->CalcLocalH();
2687 // Care of vertices internal in solids and internal faces (issue 0020676)
2688 if ( internals.hasInternalVertexInSolid() || internals.hasInternalFaces() )
2690 // store computed faces in SMESH in order not to create SMESH
2691 // faces for ng faces added here
2692 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2693 // add ng faces to solids with internal vertices
2694 AddIntVerticesInSolids( occgeo, *_ngMesh, nodeVec, internals );
2695 // duplicate mesh faces on internal faces
2696 FixIntFaces( occgeo, *_ngMesh, internals );
2697 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2699 // Let netgen compute 3D mesh
2700 startWith = endWith = netgen::MESHCONST_MESHVOLUME;
2705 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2707 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2709 if(netgen::multithread.terminate)
2712 if ( comment.empty() ) // do not overwrite a previos error
2713 comment << text(err);
2715 catch (Standard_Failure& ex)
2717 if ( comment.empty() ) // do not overwrite a previos error
2718 comment << text(ex);
2721 catch (netgen::NgException exc)
2723 if ( comment.empty() ) // do not overwrite a previos error
2724 comment << text(exc);
2727 _ticTime = ( doneTime += voluMeshingTime ) / _totalTime / _progressTic;
2729 // Let netgen optimize 3D mesh
2730 if ( !err && _optimize )
2732 startWith = endWith = netgen::MESHCONST_OPTVOLUME;
2737 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2739 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2741 if(netgen::multithread.terminate)
2744 if ( comment.empty() ) // do not overwrite a previos error
2745 comment << text(err);
2747 catch (Standard_Failure& ex)
2749 if ( comment.empty() ) // do not overwrite a previos error
2750 comment << text(ex);
2752 catch (netgen::NgException exc)
2754 if ( comment.empty() ) // do not overwrite a previos error
2755 comment << text(exc);
2759 if (!err && mparams.secondorder > 0)
2764 if ( !meshedSM[ MeshDim_1D ].empty() )
2766 // remove segments not attached to geometry (IPAL0052479)
2767 for (int i = 1; i <= _ngMesh->GetNSeg(); ++i)
2769 const netgen::Segment & seg = _ngMesh->LineSegment (i);
2770 if ( seg.epgeominfo[ 0 ].edgenr == 0 )
2771 _ngMesh->DeleteSegment( i );
2773 _ngMesh->Compress();
2775 // convert to quadratic
2776 netgen::OCCRefinementSurfaces ref (occgeo);
2777 ref.MakeSecondOrder (*_ngMesh);
2779 // care of elements already loaded to SMESH
2780 // if ( initState._nbSegments > 0 )
2781 // makeQuadratic( occgeo.emap, _mesh );
2782 // if ( initState._nbFaces > 0 )
2783 // makeQuadratic( occgeo.fmap, _mesh );
2785 catch (Standard_Failure& ex)
2787 if ( comment.empty() ) // do not overwrite a previos error
2788 comment << "Exception in netgen at passing to 2nd order ";
2790 catch (netgen::NgException exc)
2792 if ( comment.empty() ) // do not overwrite a previos error
2793 comment << exc.What();
2798 _ticTime = 0.98 / _progressTic;
2800 int nbNod = _ngMesh->GetNP();
2801 int nbSeg = _ngMesh->GetNSeg();
2802 int nbFac = _ngMesh->GetNSE();
2803 int nbVol = _ngMesh->GetNE();
2804 bool isOK = ( !err && (_isVolume ? (nbVol > 0) : (nbFac > 0)) );
2806 MESSAGE((err ? "Mesh Generation failure" : "End of Mesh Generation") <<
2807 ", nb nodes: " << nbNod <<
2808 ", nb segments: " << nbSeg <<
2809 ", nb faces: " << nbFac <<
2810 ", nb volumes: " << nbVol);
2812 // Feed back the SMESHDS with the generated Nodes and Elements
2813 if ( true /*isOK*/ ) // get whatever built
2815 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2817 if ( quadHelper.GetIsQuadratic() ) // remove free nodes
2818 for ( size_t i = 0; i < nodeVec.size(); ++i )
2819 if ( nodeVec[i] && nodeVec[i]->NbInverseElements() == 0 )
2820 _mesh->GetMeshDS()->RemoveFreeNode( nodeVec[i], 0, /*fromGroups=*/false );
2822 SMESH_ComputeErrorPtr readErr = ReadErrors(nodeVec);
2823 if ( readErr && !readErr->myBadElements.empty() )
2826 if ( !comment.empty() && !readErr->myComment.empty() ) comment += "\n";
2827 comment += readErr->myComment;
2829 if ( error->IsOK() && ( !isOK || comment.size() > 0 ))
2830 error->myName = COMPERR_ALGO_FAILED;
2831 if ( !comment.empty() )
2832 error->myComment = comment;
2834 // SetIsAlwaysComputed( true ) to empty sub-meshes, which
2835 // appear if the geometry contains coincident sub-shape due
2836 // to bool merge_solids = 1; in netgen/libsrc/occ/occgenmesh.cpp
2837 const int nbMaps = 2;
2838 const TopTools_IndexedMapOfShape* geoMaps[nbMaps] =
2839 { & occgeo.vmap, & occgeo.emap/*, & occgeo.fmap*/ };
2840 for ( int iMap = 0; iMap < nbMaps; ++iMap )
2841 for (int i = 1; i <= geoMaps[iMap]->Extent(); i++)
2842 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( geoMaps[iMap]->FindKey(i)))
2843 if ( !sm->IsMeshComputed() )
2844 sm->SetIsAlwaysComputed( true );
2846 // set bad compute error to subshapes of all failed sub-shapes
2847 if ( !error->IsOK() )
2849 bool pb2D = false, pb3D = false;
2850 for (int i = 1; i <= occgeo.fmap.Extent(); i++) {
2851 int status = occgeo.facemeshstatus[i-1];
2852 if (status == 1 ) continue;
2853 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.fmap( i ))) {
2854 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
2855 if ( !smError || smError->IsOK() ) {
2857 smError.reset( new SMESH_ComputeError( *error ));
2859 smError.reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED, "Ignored" ));
2860 if ( SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
2861 smError->myName = COMPERR_WARNING;
2863 pb2D = pb2D || smError->IsKO();
2866 if ( !pb2D ) // all faces are OK
2867 for (int i = 1; i <= occgeo.somap.Extent(); i++)
2868 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.somap( i )))
2870 bool smComputed = nbVol && !sm->IsEmpty();
2871 if ( smComputed && internals.hasInternalVertexInSolid( sm->GetId() ))
2873 int nbIntV = internals.getSolidsWithVertices().find( sm->GetId() )->second.size();
2874 SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
2875 smComputed = ( smDS->NbElements() > 0 || smDS->NbNodes() > nbIntV );
2877 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
2878 if ( !smComputed && ( !smError || smError->IsOK() ))
2880 smError.reset( new SMESH_ComputeError( *error ));
2881 if ( nbVol && SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
2883 smError->myName = COMPERR_WARNING;
2885 else if ( !smError->myBadElements.empty() ) // bad surface mesh
2887 if ( !hasBadElemOnSolid( smError->myBadElements, sm ))
2891 pb3D = pb3D || ( smError && smError->IsKO() );
2893 if ( !pb2D && !pb3D )
2894 err = 0; // no fatal errors, only warnings
2897 ngLib._isComputeOk = !err;
2902 //=============================================================================
2906 //=============================================================================
2907 bool NETGENPlugin_Mesher::Evaluate(MapShapeNbElems& aResMap)
2909 netgen::MeshingParameters& mparams = netgen::mparam;
2912 // -------------------------
2913 // Prepare OCC geometry
2914 // -------------------------
2915 netgen::OCCGeometry occgeo;
2916 list< SMESH_subMesh* > meshedSM[4]; // for 0-3 dimensions
2917 NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
2918 PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
2920 bool tooManyElems = false;
2921 const int hugeNb = std::numeric_limits<int>::max() / 100;
2926 // pass 1D simple parameters to NETGEN
2929 // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
2930 mparams.uselocalh = false;
2931 mparams.grading = 0.8; // not limitited size growth
2933 if ( _simpleHyp->GetNumberOfSegments() )
2935 mparams.maxh = occgeo.boundingbox.Diam();
2938 mparams.maxh = _simpleHyp->GetLocalLength();
2941 if ( mparams.maxh == 0.0 )
2942 mparams.maxh = occgeo.boundingbox.Diam();
2943 if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
2944 mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
2946 // let netgen create _ngMesh and calculate element size on not meshed shapes
2947 NETGENPlugin_NetgenLibWrapper ngLib;
2948 netgen::Mesh *ngMesh = NULL;
2952 int startWith = netgen::MESHCONST_ANALYSE;
2953 int endWith = netgen::MESHCONST_MESHEDGES;
2955 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, mparams, startWith, endWith);
2957 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, startWith, endWith, optstr);
2960 if(netgen::multithread.terminate)
2963 ngLib.setMesh(( Ng_Mesh*) ngMesh );
2965 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( _shape ))
2966 sm->GetComputeError().reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED ));
2971 // Pass 1D simple parameters to NETGEN
2972 // --------------------------------
2973 int nbSeg = _simpleHyp->GetNumberOfSegments();
2974 double segSize = _simpleHyp->GetLocalLength();
2975 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
2977 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
2979 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
2980 setLocalSize( e, segSize, *ngMesh );
2983 else // if ( ! _simpleHyp )
2985 // Local size on vertices and edges
2986 // --------------------------------
2987 for(std::map<int,double>::const_iterator it=EdgeId2LocalSize.begin(); it!=EdgeId2LocalSize.end(); it++)
2989 int key = (*it).first;
2990 double hi = (*it).second;
2991 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2992 const TopoDS_Edge& e = TopoDS::Edge(shape);
2993 setLocalSize( e, hi, *ngMesh );
2995 for(std::map<int,double>::const_iterator it=VertexId2LocalSize.begin(); it!=VertexId2LocalSize.end(); it++)
2997 int key = (*it).first;
2998 double hi = (*it).second;
2999 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3000 const TopoDS_Vertex& v = TopoDS::Vertex(shape);
3001 gp_Pnt p = BRep_Tool::Pnt(v);
3002 NETGENPlugin_Mesher::RestrictLocalSize( *ngMesh, p.XYZ(), hi );
3004 for(map<int,double>::const_iterator it=FaceId2LocalSize.begin();
3005 it!=FaceId2LocalSize.end(); it++)
3007 int key = (*it).first;
3008 double val = (*it).second;
3009 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3010 int faceNgID = occgeo.fmap.FindIndex(shape);
3011 occgeo.SetFaceMaxH(faceNgID, val);
3012 for ( TopExp_Explorer edgeExp( shape, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
3013 setLocalSize( TopoDS::Edge( edgeExp.Current() ), val, *ngMesh );
3016 // calculate total nb of segments and length of edges
3017 double fullLen = 0.0;
3019 int entity = mparams.secondorder > 0 ? SMDSEntity_Quad_Edge : SMDSEntity_Edge;
3020 TopTools_DataMapOfShapeInteger Edge2NbSeg;
3021 for (TopExp_Explorer exp(_shape, TopAbs_EDGE); exp.More(); exp.Next())
3023 TopoDS_Edge E = TopoDS::Edge( exp.Current() );
3024 if( !Edge2NbSeg.Bind(E,0) )
3027 double aLen = SMESH_Algo::EdgeLength(E);
3030 vector<int>& aVec = aResMap[_mesh->GetSubMesh(E)];
3032 aVec.resize( SMDSEntity_Last, 0);
3034 fullNbSeg += aVec[ entity ];
3037 // store nb of segments computed by Netgen
3038 NCollection_Map<Link> linkMap;
3039 for (int i = 1; i <= ngMesh->GetNSeg(); ++i )
3041 const netgen::Segment& seg = ngMesh->LineSegment(i);
3042 Link link(seg[0], seg[1]);
3043 if ( !linkMap.Add( link )) continue;
3044 int aGeomEdgeInd = seg.epgeominfo[0].edgenr;
3045 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
3047 vector<int>& aVec = aResMap[_mesh->GetSubMesh(occgeo.emap(aGeomEdgeInd))];
3051 // store nb of nodes on edges computed by Netgen
3052 TopTools_DataMapIteratorOfDataMapOfShapeInteger Edge2NbSegIt(Edge2NbSeg);
3053 for (; Edge2NbSegIt.More(); Edge2NbSegIt.Next())
3055 vector<int>& aVec = aResMap[_mesh->GetSubMesh(Edge2NbSegIt.Key())];
3056 if ( aVec[ entity ] > 1 && aVec[ SMDSEntity_Node ] == 0 )
3057 aVec[SMDSEntity_Node] = mparams.secondorder > 0 ? 2*aVec[ entity ]-1 : aVec[ entity ]-1;
3059 fullNbSeg += aVec[ entity ];
3060 Edge2NbSeg( Edge2NbSegIt.Key() ) = aVec[ entity ];
3062 if ( fullNbSeg == 0 )
3069 if ( double area = _simpleHyp->GetMaxElementArea() ) {
3071 mparams.maxh = sqrt(2. * area/sqrt(3.0));
3072 mparams.grading = 0.4; // moderate size growth
3075 // length from edges
3076 mparams.maxh = fullLen/fullNbSeg;
3077 mparams.grading = 0.2; // slow size growth
3080 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3081 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3083 for (TopExp_Explorer exp(_shape, TopAbs_FACE); exp.More(); exp.Next())
3085 TopoDS_Face F = TopoDS::Face( exp.Current() );
3086 SMESH_subMesh *sm = _mesh->GetSubMesh(F);
3088 BRepGProp::SurfaceProperties(F,G);
3089 double anArea = G.Mass();
3090 tooManyElems = tooManyElems || ( anArea/hugeNb > mparams.maxh*mparams.maxh );
3092 if ( !tooManyElems )
3094 TopTools_MapOfShape egdes;
3095 for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next())
3096 if ( egdes.Add( exp1.Current() ))
3097 nb1d += Edge2NbSeg.Find(exp1.Current());
3099 int nbFaces = tooManyElems ? hugeNb : int( 4*anArea / (mparams.maxh*mparams.maxh*sqrt(3.)));
3100 int nbNodes = tooManyElems ? hugeNb : (( nbFaces*3 - (nb1d-1)*2 ) / 6 + 1 );
3102 vector<int> aVec(SMDSEntity_Last, 0);
3103 if( mparams.secondorder > 0 ) {
3104 int nb1d_in = (nbFaces*3 - nb1d) / 2;
3105 aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
3106 aVec[SMDSEntity_Quad_Triangle] = nbFaces;
3109 aVec[SMDSEntity_Node] = Max ( nbNodes, 0 );
3110 aVec[SMDSEntity_Triangle] = nbFaces;
3112 aResMap[sm].swap(aVec);
3119 // pass 3D simple parameters to NETGEN
3120 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
3121 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
3123 if ( double vol = simple3d->GetMaxElementVolume() ) {
3125 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
3126 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3129 // using previous length from faces
3131 mparams.grading = 0.4;
3132 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3135 BRepGProp::VolumeProperties(_shape,G);
3136 double aVolume = G.Mass();
3137 double tetrVol = 0.1179*mparams.maxh*mparams.maxh*mparams.maxh;
3138 tooManyElems = tooManyElems || ( aVolume/hugeNb > tetrVol );
3139 int nbVols = tooManyElems ? hugeNb : int(aVolume/tetrVol);
3140 int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
3141 vector<int> aVec(SMDSEntity_Last, 0 );
3142 if ( tooManyElems ) // avoid FPE
3144 aVec[SMDSEntity_Node] = hugeNb;
3145 aVec[ mparams.secondorder > 0 ? SMDSEntity_Quad_Tetra : SMDSEntity_Tetra] = hugeNb;
3149 if( mparams.secondorder > 0 ) {
3150 aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
3151 aVec[SMDSEntity_Quad_Tetra] = nbVols;
3154 aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
3155 aVec[SMDSEntity_Tetra] = nbVols;
3158 SMESH_subMesh *sm = _mesh->GetSubMesh(_shape);
3159 aResMap[sm].swap(aVec);
3165 double NETGENPlugin_Mesher::GetProgress(const SMESH_Algo* holder,
3166 const int * algoProgressTic,
3167 const double * algoProgress) const
3169 ((int&) _progressTic ) = *algoProgressTic + 1;
3171 if ( !_occgeom ) return 0;
3173 double progress = -1;
3176 if ( _ticTime < 0 && netgen::multithread.task[0] == 'O'/*Optimizing surface*/ )
3178 ((double&) _ticTime ) = edgeFaceMeshingTime / _totalTime / _progressTic;
3180 else if ( !_optimize /*&& _occgeom->fmap.Extent() > 1*/ )
3182 int doneShapeIndex = -1;
3183 while ( doneShapeIndex+1 < _occgeom->facemeshstatus.Size() &&
3184 _occgeom->facemeshstatus[ doneShapeIndex+1 ])
3186 if ( doneShapeIndex+1 != _curShapeIndex )
3188 ((int&) _curShapeIndex) = doneShapeIndex+1;
3189 double doneShapeRate = _curShapeIndex / double( _occgeom->fmap.Extent() );
3190 double doneTime = edgeMeshingTime + doneShapeRate * faceMeshingTime;
3191 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3192 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3193 // << " " << doneTime / _totalTime / _progressTic << endl;
3197 else if ( !_optimize && _occgeom->somap.Extent() > 1 )
3199 int curShapeIndex = _curShapeIndex;
3200 if ( _ngMesh->GetNE() > 0 )
3202 netgen::Element el = (*_ngMesh)[netgen::ElementIndex( _ngMesh->GetNE()-1 )];
3203 curShapeIndex = el.GetIndex();
3205 if ( curShapeIndex != _curShapeIndex )
3207 ((int&) _curShapeIndex) = curShapeIndex;
3208 double doneShapeRate = _curShapeIndex / double( _occgeom->somap.Extent() );
3209 double doneTime = edgeFaceMeshingTime + doneShapeRate * voluMeshingTime;
3210 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3211 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3212 // << " " << doneTime / _totalTime / _progressTic << endl;
3216 progress = Max( *algoProgressTic * _ticTime, *algoProgress );
3219 ((int&) *algoProgressTic )++;
3220 ((double&) *algoProgress) = progress;
3222 //cout << progress << " " << *algoProgressTic << " " << netgen::multithread.task << " "<< _ticTime << endl;
3224 return Min( progress, 0.99 );
3227 //================================================================================
3229 * \brief Remove "test.out" and "problemfaces" files in current directory
3231 //================================================================================
3233 void NETGENPlugin_Mesher::RemoveTmpFiles()
3235 bool rm = SMESH_File("test.out").remove() ;
3237 if (rm && netgen::testout)
3239 delete netgen::testout;
3240 netgen::testout = 0;
3243 SMESH_File("problemfaces").remove();
3244 SMESH_File("occmesh.rep").remove();
3247 //================================================================================
3249 * \brief Read mesh entities preventing successful computation from "test.out" file
3251 //================================================================================
3253 SMESH_ComputeErrorPtr
3254 NETGENPlugin_Mesher::ReadErrors(const vector<const SMDS_MeshNode* >& nodeVec)
3256 SMESH_ComputeErrorPtr err = SMESH_ComputeError::New
3257 (COMPERR_BAD_INPUT_MESH, "Some edges multiple times in surface mesh");
3258 SMESH_File file("test.out");
3260 vector<int> three1(3), three2(3);
3261 const char* badEdgeStr = " multiple times in surface mesh";
3262 const int badEdgeStrLen = strlen( badEdgeStr );
3264 while( !file.eof() )
3266 if ( strncmp( file, "Edge ", 5 ) == 0 &&
3267 file.getInts( two ) &&
3268 strncmp( file, badEdgeStr, badEdgeStrLen ) == 0 &&
3269 two[0] < nodeVec.size() && two[1] < nodeVec.size())
3271 err->myBadElements.push_back( new SMDS_LinearEdge( nodeVec[ two[0]], nodeVec[ two[1]] ));
3272 file += badEdgeStrLen;
3274 else if ( strncmp( file, "Intersecting: ", 14 ) == 0 )
3277 // openelement 18 with open element 126
3281 const char* pos = file;
3282 bool ok = ( strncmp( file, "openelement ", 12 ) == 0 );
3283 ok = ok && file.getInts( two );
3284 ok = ok && file.getInts( three1 );
3285 ok = ok && file.getInts( three2 );
3286 for ( int i = 0; ok && i < 3; ++i )
3287 ok = ( three1[i] < nodeVec.size() && nodeVec[ three1[i]]);
3288 for ( int i = 0; ok && i < 3; ++i )
3289 ok = ( three2[i] < nodeVec.size() && nodeVec[ three2[i]]);
3292 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three1[0]],
3293 nodeVec[ three1[1]],
3294 nodeVec[ three1[2]]));
3295 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three2[0]],
3296 nodeVec[ three2[1]],
3297 nodeVec[ three2[2]]));
3298 err->myComment = "Intersecting triangles";
3312 size_t nbBadElems = err->myBadElements.size();
3319 //================================================================================
3321 * \brief Write a python script creating an equivalent SALOME mesh.
3322 * This is useful to see what mesh is passed as input for the next step of mesh
3323 * generation (of mesh of higher dimension)
3325 //================================================================================
3327 void NETGENPlugin_Mesher::toPython( const netgen::Mesh* ngMesh,
3328 const std::string& pyFile)
3330 ofstream outfile(pyFile.c_str(), ios::out);
3331 if ( !outfile ) return;
3333 outfile << "import SMESH" << endl
3334 << "from salome.smesh import smeshBuilder" << endl
3335 << "smesh = smeshBuilder.New(salome.myStudy)" << endl
3336 << "mesh = smesh.Mesh()" << endl << endl;
3338 using namespace netgen;
3340 for (pi = PointIndex::BASE;
3341 pi < ngMesh->GetNP()+PointIndex::BASE; pi++)
3343 outfile << "mesh.AddNode( ";
3344 outfile << (*ngMesh)[pi](0) << ", ";
3345 outfile << (*ngMesh)[pi](1) << ", ";
3346 outfile << (*ngMesh)[pi](2) << ") ## "<< pi << endl;
3349 int nbDom = ngMesh->GetNDomains();
3350 for ( int i = 0; i < nbDom; ++i )
3351 outfile<< "grp" << i+1 << " = mesh.CreateEmptyGroup( SMESH.FACE, 'domain"<< i+1 << "')"<< endl;
3353 SurfaceElementIndex sei;
3354 for (sei = 0; sei < ngMesh->GetNSE(); sei++)
3356 outfile << "mesh.AddFace([ ";
3357 Element2d sel = (*ngMesh)[sei];
3358 for (int j = 0; j < sel.GetNP(); j++)
3359 outfile << sel[j] << ( j+1 < sel.GetNP() ? ", " : " ])");
3360 if ( sel.IsDeleted() ) outfile << " ## IsDeleted ";
3363 if ((*ngMesh)[sei].GetIndex())
3365 if ( int dom1 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainIn())
3366 outfile << "grp"<< dom1 <<".Add([ " << (int)sei+1 << " ])" << endl;
3367 if ( int dom2 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainOut())
3368 outfile << "grp"<< dom2 <<".Add([ " << (int)sei+1 << " ])" << endl;
3372 for (ElementIndex ei = 0; ei < ngMesh->GetNE(); ei++)
3374 Element el = (*ngMesh)[ei];
3375 outfile << "mesh.AddVolume([ ";
3376 for (int j = 0; j < el.GetNP(); j++)
3377 outfile << el[j] << ( j+1 < el.GetNP() ? ", " : " ])");
3381 for (int i = 1; i <= ngMesh->GetNSeg(); i++)
3383 const Segment & seg = ngMesh->LineSegment (i);
3384 outfile << "mesh.AddEdge([ "
3386 << seg[1] << " ])" << endl;
3388 cout << "Write " << pyFile << endl;
3391 //================================================================================
3393 * \brief Constructor of NETGENPlugin_ngMeshInfo
3395 //================================================================================
3397 NETGENPlugin_ngMeshInfo::NETGENPlugin_ngMeshInfo( netgen::Mesh* ngMesh):
3402 _nbNodes = ngMesh->GetNP();
3403 _nbSegments = ngMesh->GetNSeg();
3404 _nbFaces = ngMesh->GetNSE();
3405 _nbVolumes = ngMesh->GetNE();
3409 _nbNodes = _nbSegments = _nbFaces = _nbVolumes = 0;
3413 //================================================================================
3415 * \brief Copy LocalH member from one netgen mesh to another
3417 //================================================================================
3419 void NETGENPlugin_ngMeshInfo::transferLocalH( netgen::Mesh* fromMesh,
3420 netgen::Mesh* toMesh )
3422 if ( !fromMesh->LocalHFunctionGenerated() ) return;
3423 if ( !toMesh->LocalHFunctionGenerated() )
3425 toMesh->CalcLocalH(netgen::mparam.grading);
3427 toMesh->CalcLocalH();
3430 const size_t size = sizeof( netgen::LocalH );
3431 _copyOfLocalH = new char[ size ];
3432 memcpy( (void*)_copyOfLocalH, (void*)&toMesh->LocalHFunction(), size );
3433 memcpy( (void*)&toMesh->LocalHFunction(), (void*)&fromMesh->LocalHFunction(), size );
3436 //================================================================================
3438 * \brief Restore LocalH member of a netgen mesh
3440 //================================================================================
3442 void NETGENPlugin_ngMeshInfo::restoreLocalH( netgen::Mesh* toMesh )
3444 if ( _copyOfLocalH )
3446 const size_t size = sizeof( netgen::LocalH );
3447 memcpy( (void*)&toMesh->LocalHFunction(), (void*)_copyOfLocalH, size );
3448 delete [] _copyOfLocalH;
3453 //================================================================================
3455 * \brief Find "internal" sub-shapes
3457 //================================================================================
3459 NETGENPlugin_Internals::NETGENPlugin_Internals( SMESH_Mesh& mesh,
3460 const TopoDS_Shape& shape,
3462 : _mesh( mesh ), _is3D( is3D )
3464 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3466 TopExp_Explorer f,e;
3467 for ( f.Init( shape, TopAbs_FACE ); f.More(); f.Next() )
3469 int faceID = meshDS->ShapeToIndex( f.Current() );
3471 // find not computed internal edges
3473 for ( e.Init( f.Current().Oriented(TopAbs_FORWARD), TopAbs_EDGE ); e.More(); e.Next() )
3474 if ( e.Current().Orientation() == TopAbs_INTERNAL )
3476 SMESH_subMesh* eSM = mesh.GetSubMesh( e.Current() );
3477 if ( eSM->IsEmpty() )
3479 _e2face.insert( make_pair( eSM->GetId(), faceID ));
3480 for ( TopoDS_Iterator v(e.Current()); v.More(); v.Next() )
3481 _e2face.insert( make_pair( meshDS->ShapeToIndex( v.Value() ), faceID ));
3485 // find internal vertices in a face
3486 set<int> intVV; // issue 0020850 where same vertex is twice in a face
3487 for ( TopoDS_Iterator fSub( f.Current() ); fSub.More(); fSub.Next())
3488 if ( fSub.Value().ShapeType() == TopAbs_VERTEX )
3490 int vID = meshDS->ShapeToIndex( fSub.Value() );
3491 if ( intVV.insert( vID ).second )
3492 _f2v[ faceID ].push_back( vID );
3497 // find internal faces and their subshapes where nodes are to be doubled
3498 // to make a crack with non-sewed borders
3500 if ( f.Current().Orientation() == TopAbs_INTERNAL )
3502 _intShapes.insert( meshDS->ShapeToIndex( f.Current() ));
3505 list< TopoDS_Shape > edges;
3506 for ( e.Init( f.Current(), TopAbs_EDGE ); e.More(); e.Next())
3507 if ( SMESH_MesherHelper::NbAncestors( e.Current(), mesh, TopAbs_FACE ) > 1 )
3509 _intShapes.insert( meshDS->ShapeToIndex( e.Current() ));
3510 edges.push_back( e.Current() );
3511 // find border faces
3512 PShapeIteratorPtr fIt =
3513 SMESH_MesherHelper::GetAncestors( edges.back(),mesh,TopAbs_FACE );
3514 while ( const TopoDS_Shape* pFace = fIt->next() )
3515 if ( !pFace->IsSame( f.Current() ))
3516 _borderFaces.insert( meshDS->ShapeToIndex( *pFace ));
3519 // we consider vertex internal if it is shared by more than one internal edge
3520 list< TopoDS_Shape >::iterator edge = edges.begin();
3521 for ( ; edge != edges.end(); ++edge )
3522 for ( TopoDS_Iterator v( *edge ); v.More(); v.Next() )
3524 set<int> internalEdges;
3525 PShapeIteratorPtr eIt =
3526 SMESH_MesherHelper::GetAncestors( v.Value(),mesh,TopAbs_EDGE );
3527 while ( const TopoDS_Shape* pEdge = eIt->next() )
3529 int edgeID = meshDS->ShapeToIndex( *pEdge );
3530 if ( isInternalShape( edgeID ))
3531 internalEdges.insert( edgeID );
3533 if ( internalEdges.size() > 1 )
3534 _intShapes.insert( meshDS->ShapeToIndex( v.Value() ));
3538 } // loop on geom faces
3540 // find vertices internal in solids
3543 for ( TopExp_Explorer so(shape, TopAbs_SOLID); so.More(); so.Next())
3545 int soID = meshDS->ShapeToIndex( so.Current() );
3546 for ( TopoDS_Iterator soSub( so.Current() ); soSub.More(); soSub.Next())
3547 if ( soSub.Value().ShapeType() == TopAbs_VERTEX )
3548 _s2v[ soID ].push_back( meshDS->ShapeToIndex( soSub.Value() ));
3553 //================================================================================
3555 * \brief Find mesh faces on non-internal geom faces sharing internal edge
3556 * some nodes of which are to be doubled to make the second border of the "crack"
3558 //================================================================================
3560 void NETGENPlugin_Internals::findBorderElements( TIDSortedElemSet & borderElems )
3562 if ( _intShapes.empty() ) return;
3564 SMESH_Mesh& mesh = const_cast<SMESH_Mesh&>(_mesh);
3565 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3567 // loop on internal geom edges
3568 set<int>::const_iterator intShapeId = _intShapes.begin();
3569 for ( ; intShapeId != _intShapes.end(); ++intShapeId )
3571 const TopoDS_Shape& s = meshDS->IndexToShape( *intShapeId );
3572 if ( s.ShapeType() != TopAbs_EDGE ) continue;
3574 // get internal and non-internal geom faces sharing the internal edge <s>
3576 set<int>::iterator bordFace = _borderFaces.end();
3577 PShapeIteratorPtr faces = SMESH_MesherHelper::GetAncestors( s, _mesh, TopAbs_FACE );
3578 while ( const TopoDS_Shape* pFace = faces->next() )
3580 int faceID = meshDS->ShapeToIndex( *pFace );
3581 if ( isInternalShape( faceID ))
3584 bordFace = _borderFaces.insert( faceID ).first;
3586 if ( bordFace == _borderFaces.end() || !intFace ) continue;
3588 // get all links of mesh faces on internal geom face sharing nodes on edge <s>
3589 set< SMESH_OrientedLink > links; //!< links of faces on internal geom face
3590 list<const SMDS_MeshElement*> suspectFaces[2]; //!< mesh faces on border geom faces
3591 int nbSuspectFaces = 0;
3592 SMESHDS_SubMesh* intFaceSM = meshDS->MeshElements( intFace );
3593 if ( !intFaceSM || intFaceSM->NbElements() == 0 ) continue;
3594 SMESH_subMeshIteratorPtr smIt = mesh.GetSubMesh( s )->getDependsOnIterator(true,true);
3595 while ( smIt->more() )
3597 SMESHDS_SubMesh* sm = smIt->next()->GetSubMeshDS();
3598 if ( !sm ) continue;
3599 SMDS_NodeIteratorPtr nIt = sm->GetNodes();
3600 while ( nIt->more() )
3602 const SMDS_MeshNode* nOnEdge = nIt->next();
3603 SMDS_ElemIteratorPtr fIt = nOnEdge->GetInverseElementIterator(SMDSAbs_Face);
3604 while ( fIt->more() )
3606 const SMDS_MeshElement* f = fIt->next();
3607 int nbNodes = f->NbNodes() / ( f->IsQuadratic() ? 2 : 1 );
3608 if ( intFaceSM->Contains( f ))
3610 for ( int i = 0; i < nbNodes; ++i )
3611 links.insert( SMESH_OrientedLink( f->GetNode(i), f->GetNode((i+1)%nbNodes)));
3616 for ( int i = 0; i < nbNodes; ++i )
3617 nbDblNodes += isInternalShape( f->GetNode(i)->getshapeId() );
3619 suspectFaces[ nbDblNodes < 2 ].push_back( f );
3625 // suspectFaces[0] having link with same orientation as mesh faces on
3626 // the internal geom face are <borderElems>. suspectFaces[1] have
3627 // only one node on edge <s>, we decide on them later (at the 2nd loop)
3628 // by links of <borderElems> found at the 1st and 2nd loops
3629 set< SMESH_OrientedLink > borderLinks;
3630 for ( int isPostponed = 0; isPostponed < 2; ++isPostponed )
3632 list<const SMDS_MeshElement*>::iterator fIt = suspectFaces[isPostponed].begin();
3633 for ( int nbF = 0; fIt != suspectFaces[isPostponed].end(); ++fIt, ++nbF )
3635 const SMDS_MeshElement* f = *fIt;
3636 bool isBorder = false, linkFound = false, borderLinkFound = false;
3637 list< SMESH_OrientedLink > faceLinks;
3638 int nbNodes = f->NbNodes() / ( f->IsQuadratic() ? 2 : 1 );
3639 for ( int i = 0; i < nbNodes; ++i )
3641 SMESH_OrientedLink link( f->GetNode(i), f->GetNode((i+1)%nbNodes));
3642 faceLinks.push_back( link );
3645 set< SMESH_OrientedLink >::iterator foundLink = links.find( link );
3646 if ( foundLink != links.end() )
3649 isBorder = ( foundLink->_reversed == link._reversed );
3650 if ( !isBorder && !isPostponed ) break;
3651 faceLinks.pop_back();
3653 else if ( isPostponed && !borderLinkFound )
3655 foundLink = borderLinks.find( link );
3656 if ( foundLink != borderLinks.end() )
3658 borderLinkFound = true;
3659 isBorder = ( foundLink->_reversed != link._reversed );
3666 borderElems.insert( f );
3667 borderLinks.insert( faceLinks.begin(), faceLinks.end() );
3669 else if ( !linkFound && !borderLinkFound )
3671 suspectFaces[1].push_back( f );
3672 if ( nbF > 2 * nbSuspectFaces )
3673 break; // dead loop protection
3680 //================================================================================
3682 * \brief put internal shapes in maps and fill in submeshes to precompute
3684 //================================================================================
3686 void NETGENPlugin_Internals::getInternalEdges( TopTools_IndexedMapOfShape& fmap,
3687 TopTools_IndexedMapOfShape& emap,
3688 TopTools_IndexedMapOfShape& vmap,
3689 list< SMESH_subMesh* > smToPrecompute[])
3691 if ( !hasInternalEdges() ) return;
3692 map<int,int>::const_iterator ev_face = _e2face.begin();
3693 for ( ; ev_face != _e2face.end(); ++ev_face )
3695 const TopoDS_Shape& ev = _mesh.GetMeshDS()->IndexToShape( ev_face->first );
3696 const TopoDS_Shape& face = _mesh.GetMeshDS()->IndexToShape( ev_face->second );
3698 ( ev.ShapeType() == TopAbs_EDGE ? emap : vmap ).Add( ev );
3700 //cout<<"INTERNAL EDGE or VERTEX "<<ev_face->first<<" on face "<<ev_face->second<<endl;
3702 smToPrecompute[ MeshDim_1D ].push_back( _mesh.GetSubMeshContaining( ev_face->first ));
3706 //================================================================================
3708 * \brief return shapes and submeshes to be meshed and already meshed boundary submeshes
3710 //================================================================================
3712 void NETGENPlugin_Internals::getInternalFaces( TopTools_IndexedMapOfShape& fmap,
3713 TopTools_IndexedMapOfShape& emap,
3714 list< SMESH_subMesh* >& intFaceSM,
3715 list< SMESH_subMesh* >& boundarySM)
3717 if ( !hasInternalFaces() ) return;
3719 // <fmap> and <emap> are for not yet meshed shapes
3720 // <intFaceSM> is for submeshes of faces
3721 // <boundarySM> is for meshed edges and vertices
3726 set<int> shapeIDs ( _intShapes );
3727 if ( !_borderFaces.empty() )
3728 shapeIDs.insert( _borderFaces.begin(), _borderFaces.end() );
3730 set<int>::const_iterator intS = shapeIDs.begin();
3731 for ( ; intS != shapeIDs.end(); ++intS )
3733 SMESH_subMesh* sm = _mesh.GetSubMeshContaining( *intS );
3735 if ( sm->GetSubShape().ShapeType() != TopAbs_FACE ) continue;
3737 intFaceSM.push_back( sm );
3739 // add submeshes of not computed internal faces
3740 if ( !sm->IsEmpty() ) continue;
3742 SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(true,true);
3743 while ( smIt->more() )
3746 const TopoDS_Shape& s = sm->GetSubShape();
3748 if ( sm->IsEmpty() )
3751 switch ( s.ShapeType() ) {
3752 case TopAbs_FACE: fmap.Add ( s ); break;
3753 case TopAbs_EDGE: emap.Add ( s ); break;
3759 if ( s.ShapeType() != TopAbs_FACE )
3760 boundarySM.push_back( sm );
3766 //================================================================================
3768 * \brief Return true if given shape is to be precomputed in order to be correctly
3769 * added to netgen mesh
3771 //================================================================================
3773 bool NETGENPlugin_Internals::isShapeToPrecompute(const TopoDS_Shape& s)
3775 int shapeID = _mesh.GetMeshDS()->ShapeToIndex( s );
3776 switch ( s.ShapeType() ) {
3777 case TopAbs_FACE : break; //return isInternalShape( shapeID ) || isBorderFace( shapeID );
3778 case TopAbs_EDGE : return isInternalEdge( shapeID );
3779 case TopAbs_VERTEX: break;
3785 //================================================================================
3787 * \brief Return SMESH
3789 //================================================================================
3791 SMESH_Mesh& NETGENPlugin_Internals::getMesh() const
3793 return const_cast<SMESH_Mesh&>( _mesh );
3796 //================================================================================
3798 * \brief Initialize netgen library
3800 //================================================================================
3802 NETGENPlugin_NetgenLibWrapper::NETGENPlugin_NetgenLibWrapper()
3806 _isComputeOk = false;
3808 if ( !getenv( "KEEP_NETGEN_OUTPUT" ))
3810 // redirect all netgen output (mycout,myerr,cout) to _outputFileName
3811 _outputFileName = getOutputFileName();
3812 netgen::mycout = new ofstream ( _outputFileName.c_str() );
3813 netgen::myerr = netgen::mycout;
3814 _coutBuffer = std::cout.rdbuf();
3816 cout << "NOTE: netgen output is redirected to file " << _outputFileName << endl;
3818 std::cout.rdbuf( netgen::mycout->rdbuf() );
3822 _ngMesh = Ng_NewMesh();
3825 //================================================================================
3827 * \brief Finish using netgen library
3829 //================================================================================
3831 NETGENPlugin_NetgenLibWrapper::~NETGENPlugin_NetgenLibWrapper()
3833 Ng_DeleteMesh( _ngMesh );
3835 NETGENPlugin_Mesher::RemoveTmpFiles();
3837 std::cout.rdbuf( _coutBuffer );
3844 //================================================================================
3846 * \brief Set netgen mesh to delete at destruction
3848 //================================================================================
3850 void NETGENPlugin_NetgenLibWrapper::setMesh( Ng_Mesh* mesh )
3853 Ng_DeleteMesh( _ngMesh );
3857 //================================================================================
3859 * \brief Return a unique file name
3861 //================================================================================
3863 std::string NETGENPlugin_NetgenLibWrapper::getOutputFileName()
3865 std::string aTmpDir = SALOMEDS_Tool::GetTmpDir();
3867 TCollection_AsciiString aGenericName = (char*)aTmpDir.c_str();
3868 aGenericName += "NETGEN_";
3870 aGenericName += getpid();
3872 aGenericName += _getpid();
3874 aGenericName += "_";
3875 aGenericName += Abs((Standard_Integer)(long) aGenericName.ToCString());
3876 aGenericName += ".out";
3878 return aGenericName.ToCString();
3881 //================================================================================
3883 * \brief Remove file with netgen output
3885 //================================================================================
3887 void NETGENPlugin_NetgenLibWrapper::removeOutputFile()
3889 if ( !_outputFileName.empty() )
3891 if ( netgen::mycout )
3893 delete netgen::mycout;
3897 string tmpDir = SALOMEDS_Tool::GetDirFromPath ( _outputFileName );
3898 string aFileName = SALOMEDS_Tool::GetNameFromPath( _outputFileName ) + ".out";
3899 SALOMEDS::ListOfFileNames_var aFiles = new SALOMEDS::ListOfFileNames;
3901 aFiles[0] = aFileName.c_str();
3903 SALOMEDS_Tool::RemoveTemporaryFiles( tmpDir.c_str(), aFiles.in(), true );