1 // Copyright (C) 2007-2016 CEA/DEN, EDF R&D, OPEN CASCADE
3 // Copyright (C) 2003-2007 OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
6 // This library is free software; you can redistribute it and/or
7 // modify it under the terms of the GNU Lesser General Public
8 // License as published by the Free Software Foundation; either
9 // version 2.1 of the License, or (at your option) any later version.
11 // This library is distributed in the hope that it will be useful,
12 // but WITHOUT ANY WARRANTY; without even the implied warranty of
13 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
14 // Lesser General Public License for more details.
16 // You should have received a copy of the GNU Lesser General Public
17 // License along with this library; if not, write to the Free Software
18 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
23 // NETGENPlugin : C++ implementation
24 // File : NETGENPlugin_Mesher.cxx
25 // Author : Michael Sazonov (OCN)
28 //=============================================================================
30 #include "NETGENPlugin_Mesher.hxx"
31 #include "NETGENPlugin_Hypothesis_2D.hxx"
32 #include "NETGENPlugin_SimpleHypothesis_3D.hxx"
34 #include <SMDS_FaceOfNodes.hxx>
35 #include <SMDS_MeshElement.hxx>
36 #include <SMDS_MeshNode.hxx>
37 #include <SMESHDS_Mesh.hxx>
38 #include <SMESH_Block.hxx>
39 #include <SMESH_Comment.hxx>
40 #include <SMESH_ComputeError.hxx>
41 #include <SMESH_File.hxx>
42 #include <SMESH_Gen_i.hxx>
43 #include <SMESH_Mesh.hxx>
44 #include <SMESH_MesherHelper.hxx>
45 #include <SMESH_subMesh.hxx>
46 #include <StdMeshers_QuadToTriaAdaptor.hxx>
47 #include <StdMeshers_ViscousLayers2D.hxx>
49 #include <SALOMEDS_Tool.hxx>
51 #include <utilities.h>
53 #include <BRepBuilderAPI_Copy.hxx>
54 #include <BRep_Tool.hxx>
55 #include <Bnd_B3d.hxx>
56 #include <NCollection_Map.hxx>
57 #include <Standard_ErrorHandler.hxx>
58 #include <Standard_ProgramError.hxx>
59 #include <TColStd_MapOfInteger.hxx>
61 #include <TopExp_Explorer.hxx>
62 #include <TopTools_DataMapIteratorOfDataMapOfShapeInteger.hxx>
63 #include <TopTools_DataMapIteratorOfDataMapOfShapeShape.hxx>
64 #include <TopTools_DataMapOfShapeInteger.hxx>
65 #include <TopTools_DataMapOfShapeShape.hxx>
66 #include <TopTools_MapOfShape.hxx>
69 // Netgen include files
73 #include <occgeom.hpp>
74 #include <meshing.hpp>
75 //#include <ngexception.hpp>
78 extern int OCCGenerateMesh (OCCGeometry&, Mesh*&, MeshingParameters&, int, int);
80 extern int OCCGenerateMesh (OCCGeometry&, Mesh*&, int, int, char*);
82 //extern void OCCSetLocalMeshSize(OCCGeometry & geom, Mesh & mesh);
83 extern MeshingParameters mparam;
84 extern volatile multithreadt multithread;
85 extern bool merge_solids;
87 // values used for occgeo.facemeshstatus
88 enum EFaceMeshStatus { FACE_NOT_TREATED = 0,
100 using namespace nglib;
104 #define nodeVec_ACCESS(index) ((SMDS_MeshNode*) nodeVec.at((index)))
106 #define nodeVec_ACCESS(index) ((SMDS_MeshNode*) nodeVec[index])
109 #define NGPOINT_COORDS(p) p(0),p(1),p(2)
112 // dump elements added to ng mesh
113 //#define DUMP_SEGMENTS
114 //#define DUMP_TRIANGLES
115 //#define DUMP_TRIANGLES_SCRIPT "/tmp/trias.py" //!< debug AddIntVerticesInSolids()
118 TopTools_IndexedMapOfShape ShapesWithLocalSize;
119 std::map<int,double> VertexId2LocalSize;
120 std::map<int,double> EdgeId2LocalSize;
121 std::map<int,double> FaceId2LocalSize;
123 //=============================================================================
127 //=============================================================================
129 NETGENPlugin_Mesher::NETGENPlugin_Mesher (SMESH_Mesh* mesh,
130 const TopoDS_Shape& aShape,
136 _fineness(NETGENPlugin_Hypothesis::GetDefaultFineness()),
137 _isViscousLayers2D(false),
146 SetDefaultParameters();
147 ShapesWithLocalSize.Clear();
148 VertexId2LocalSize.clear();
149 EdgeId2LocalSize.clear();
150 FaceId2LocalSize.clear();
153 //================================================================================
157 //================================================================================
159 NETGENPlugin_Mesher::~NETGENPlugin_Mesher()
167 //================================================================================
169 * Set pointer to NETGENPlugin_Mesher* field of the holder, that will be
170 * nullified at destruction of this
172 //================================================================================
174 void NETGENPlugin_Mesher::SetSelfPointer( NETGENPlugin_Mesher ** ptr )
185 //================================================================================
187 * \brief Initialize global NETGEN parameters with default values
189 //================================================================================
191 void NETGENPlugin_Mesher::SetDefaultParameters()
193 netgen::MeshingParameters& mparams = netgen::mparam;
194 // maximal mesh edge size
195 mparams.maxh = 0;//NETGENPlugin_Hypothesis::GetDefaultMaxSize();
197 // minimal number of segments per edge
198 mparams.segmentsperedge = NETGENPlugin_Hypothesis::GetDefaultNbSegPerEdge();
199 // rate of growth of size between elements
200 mparams.grading = NETGENPlugin_Hypothesis::GetDefaultGrowthRate();
201 // safety factor for curvatures (elements per radius)
202 mparams.curvaturesafety = NETGENPlugin_Hypothesis::GetDefaultNbSegPerRadius();
203 // create elements of second order
204 mparams.secondorder = NETGENPlugin_Hypothesis::GetDefaultSecondOrder();
205 // quad-dominated surface meshing
209 mparams.quad = NETGENPlugin_Hypothesis_2D::GetDefaultQuadAllowed();
210 _fineness = NETGENPlugin_Hypothesis::GetDefaultFineness();
211 mparams.uselocalh = NETGENPlugin_Hypothesis::GetDefaultSurfaceCurvature();
212 netgen::merge_solids = NETGENPlugin_Hypothesis::GetDefaultFuseEdges();
215 //=============================================================================
219 //=============================================================================
221 void SetLocalSize(TopoDS_Shape GeomShape, double LocalSize)
223 if ( GeomShape.IsNull() ) return;
224 TopAbs_ShapeEnum GeomType = GeomShape.ShapeType();
225 if (GeomType == TopAbs_COMPOUND) {
226 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()) {
227 SetLocalSize(it.Value(), LocalSize);
232 if (! ShapesWithLocalSize.Contains(GeomShape))
233 key = ShapesWithLocalSize.Add(GeomShape);
235 key = ShapesWithLocalSize.FindIndex(GeomShape);
236 if (GeomType == TopAbs_VERTEX) {
237 VertexId2LocalSize[key] = LocalSize;
238 } else if (GeomType == TopAbs_EDGE) {
239 EdgeId2LocalSize[key] = LocalSize;
240 } else if (GeomType == TopAbs_FACE) {
241 FaceId2LocalSize[key] = LocalSize;
245 //=============================================================================
247 * Pass parameters to NETGEN
249 //=============================================================================
250 void NETGENPlugin_Mesher::SetParameters(const NETGENPlugin_Hypothesis* hyp)
254 netgen::MeshingParameters& mparams = netgen::mparam;
255 // Initialize global NETGEN parameters:
256 // maximal mesh segment size
257 mparams.maxh = hyp->GetMaxSize();
258 // maximal mesh element linear size
259 mparams.minh = hyp->GetMinSize();
260 // minimal number of segments per edge
261 mparams.segmentsperedge = hyp->GetNbSegPerEdge();
262 // rate of growth of size between elements
263 mparams.grading = hyp->GetGrowthRate();
264 // safety factor for curvatures (elements per radius)
265 mparams.curvaturesafety = hyp->GetNbSegPerRadius();
266 // create elements of second order
267 mparams.secondorder = hyp->GetSecondOrder() ? 1 : 0;
268 // quad-dominated surface meshing
269 // only triangles are allowed for volumic mesh (before realizing IMP 0021676)
271 mparams.quad = hyp->GetQuadAllowed() ? 1 : 0;
272 _optimize = hyp->GetOptimize();
273 _fineness = hyp->GetFineness();
274 mparams.uselocalh = hyp->GetSurfaceCurvature();
275 netgen::merge_solids = hyp->GetFuseEdges();
278 SMESH_Gen_i* smeshGen_i = SMESH_Gen_i::GetSMESHGen();
279 CORBA::Object_var anObject = smeshGen_i->GetNS()->Resolve("/myStudyManager");
280 SALOMEDS::StudyManager_var aStudyMgr = SALOMEDS::StudyManager::_narrow(anObject);
281 SALOMEDS::Study_var myStudy = aStudyMgr->GetStudyByID(hyp->GetStudyId());
283 const NETGENPlugin_Hypothesis::TLocalSize localSizes = hyp->GetLocalSizesAndEntries();
284 NETGENPlugin_Hypothesis::TLocalSize::const_iterator it = localSizes.begin();
285 for ( ; it != localSizes.end() ; it++)
287 std::string entry = (*it).first;
288 double val = (*it).second;
290 GEOM::GEOM_Object_var aGeomObj;
291 TopoDS_Shape S = TopoDS_Shape();
292 SALOMEDS::SObject_var aSObj = myStudy->FindObjectID( entry.c_str() );
293 if (!aSObj->_is_nil()) {
294 CORBA::Object_var obj = aSObj->GetObject();
295 aGeomObj = GEOM::GEOM_Object::_narrow(obj);
298 if ( !aGeomObj->_is_nil() )
299 S = smeshGen_i->GeomObjectToShape( aGeomObj.in() );
301 SetLocalSize(S, val);
306 //=============================================================================
308 * Pass simple parameters to NETGEN
310 //=============================================================================
312 void NETGENPlugin_Mesher::SetParameters(const NETGENPlugin_SimpleHypothesis_2D* hyp)
316 SetDefaultParameters();
319 //=============================================================================
321 * Link - a pair of integer numbers
323 //=============================================================================
327 Link(int _n1, int _n2) : n1(_n1), n2(_n2) {}
328 Link() : n1(0), n2(0) {}
329 bool Contains( int n ) const { return n == n1 || n == n2; }
330 bool IsConnected( const Link& other ) const
332 return (( Contains( other.n1 ) || Contains( other.n2 )) && ( this != &other ));
336 int HashCode(const Link& aLink, int aLimit)
338 return HashCode(aLink.n1 + aLink.n2, aLimit);
341 Standard_Boolean IsEqual(const Link& aLink1, const Link& aLink2)
343 return (( aLink1.n1 == aLink2.n1 && aLink1.n2 == aLink2.n2 ) ||
344 ( aLink1.n1 == aLink2.n2 && aLink1.n2 == aLink2.n1 ));
349 //================================================================================
351 * \brief return id of netgen point corresponding to SMDS node
353 //================================================================================
354 typedef map< const SMDS_MeshNode*, int > TNode2IdMap;
356 int ngNodeId( const SMDS_MeshNode* node,
357 netgen::Mesh& ngMesh,
358 TNode2IdMap& nodeNgIdMap)
360 int newNgId = ngMesh.GetNP() + 1;
362 TNode2IdMap::iterator node_id = nodeNgIdMap.insert( make_pair( node, newNgId )).first;
364 if ( node_id->second == newNgId)
366 #if defined(DUMP_SEGMENTS) || defined(DUMP_TRIANGLES)
367 cout << "Ng " << newNgId << " - " << node;
369 netgen::MeshPoint p( netgen::Point<3> (node->X(), node->Y(), node->Z()) );
370 ngMesh.AddPoint( p );
372 return node_id->second;
375 //================================================================================
377 * \brief Return computed EDGEs connected to the given one
379 //================================================================================
381 list< TopoDS_Edge > getConnectedEdges( const TopoDS_Edge& edge,
382 const TopoDS_Face& face,
383 const set< SMESH_subMesh* > & computedSM,
384 const SMESH_MesherHelper& helper,
385 map< SMESH_subMesh*, set< int > >& addedEdgeSM2Faces)
388 list< TopoDS_Edge > edges;
389 list< int > nbEdgesInWire;
390 /*int nbWires =*/ SMESH_Block::GetOrderedEdges( face, edges, nbEdgesInWire);
392 // find <edge> within <edges>
393 list< TopoDS_Edge >::iterator eItFwd = edges.begin();
394 for ( ; eItFwd != edges.end(); ++eItFwd )
395 if ( edge.IsSame( *eItFwd ))
397 if ( eItFwd == edges.end()) return list< TopoDS_Edge>();
399 if ( eItFwd->Orientation() >= TopAbs_INTERNAL )
401 // connected INTERNAL edges returned from GetOrderedEdges() are wrongly oriented
402 // so treat each INTERNAL edge separately
403 TopoDS_Edge e = *eItFwd;
405 edges.push_back( e );
409 // get all computed EDGEs connected to <edge>
411 list< TopoDS_Edge >::iterator eItBack = eItFwd, ePrev;
412 TopoDS_Vertex vCommon;
413 TopTools_MapOfShape eAdded; // map used not to add a seam edge twice to <edges>
416 // put edges before <edge> to <edges> back
417 while ( edges.begin() != eItFwd )
418 edges.splice( edges.end(), edges, edges.begin() );
422 while ( ++eItFwd != edges.end() )
424 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( *eItFwd );
426 bool connected = TopExp::CommonVertex( *ePrev, *eItFwd, vCommon );
427 bool computed = sm->IsMeshComputed();
428 bool added = addedEdgeSM2Faces[ sm ].count( helper.GetSubShapeID() );
429 bool doubled = !eAdded.Add( *eItFwd );
430 bool orientOK = (( ePrev ->Orientation() < TopAbs_INTERNAL ) ==
431 ( eItFwd->Orientation() < TopAbs_INTERNAL ) );
432 if ( !connected || !computed || !orientOK || added || doubled )
434 // stop advancement; move edges from tail to head
435 while ( edges.back() != *ePrev )
436 edges.splice( edges.begin(), edges, --edges.end() );
442 while ( eItBack != edges.begin() )
446 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( *eItBack );
448 bool connected = TopExp::CommonVertex( *ePrev, *eItBack, vCommon );
449 bool computed = sm->IsMeshComputed();
450 bool added = addedEdgeSM2Faces[ sm ].count( helper.GetSubShapeID() );
451 bool doubled = !eAdded.Add( *eItBack );
452 bool orientOK = (( ePrev ->Orientation() < TopAbs_INTERNAL ) ==
453 ( eItBack->Orientation() < TopAbs_INTERNAL ) );
454 if ( !connected || !computed || !orientOK || added || doubled)
457 edges.erase( edges.begin(), ePrev );
461 if ( edges.front() != edges.back() )
463 // assure that the 1st vertex is meshed
464 TopoDS_Edge eLast = edges.back();
465 while ( !SMESH_Algo::VertexNode( SMESH_MesherHelper::IthVertex( 0, edges.front()), helper.GetMeshDS())
467 edges.front() != eLast )
468 edges.splice( edges.end(), edges, edges.begin() );
473 //================================================================================
475 * \brief Make triangulation of a shape precise enough
477 //================================================================================
479 void updateTriangulation( const TopoDS_Shape& shape )
481 // static set< Poly_Triangulation* > updated;
483 // TopLoc_Location loc;
484 // TopExp_Explorer fExp( shape, TopAbs_FACE );
485 // for ( ; fExp.More(); fExp.Next() )
487 // Handle(Poly_Triangulation) triangulation =
488 // BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
489 // if ( triangulation.IsNull() ||
490 // updated.insert( triangulation.operator->() ).second )
492 // BRepTools::Clean (shape);
495 BRepMesh_IncrementalMesh e(shape, 0.01, true);
497 catch (Standard_Failure)
500 // updated.erase( triangulation.operator->() );
501 // triangulation = BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
502 // updated.insert( triangulation.operator->() );
506 //================================================================================
508 * \brief Returns a medium node either existing in SMESH of created by NETGEN
509 * \param [in] corner1 - corner node 1
510 * \param [in] corner2 - corner node 2
511 * \param [in] defaultMedium - the node created by NETGEN
512 * \param [in] helper - holder of medium nodes existing in SMESH
513 * \return const SMDS_MeshNode* - the result node
515 //================================================================================
517 const SMDS_MeshNode* mediumNode( const SMDS_MeshNode* corner1,
518 const SMDS_MeshNode* corner2,
519 const SMDS_MeshNode* defaultMedium,
520 const SMESH_MesherHelper* helper)
524 TLinkNodeMap::const_iterator l2n =
525 helper->GetTLinkNodeMap().find( SMESH_TLink( corner1, corner2 ));
526 if ( l2n != helper->GetTLinkNodeMap().end() )
527 defaultMedium = l2n->second;
529 return defaultMedium;
532 //================================================================================
534 * \brief Assure that mesh on given shapes is quadratic
536 //================================================================================
538 void makeQuadratic( const TopTools_IndexedMapOfShape& shapes,
541 for ( int i = 1; i <= shapes.Extent(); ++i )
543 SMESHDS_SubMesh* smDS = mesh->GetMeshDS()->MeshElements( shapes(i) );
544 if ( !smDS ) continue;
545 SMDS_ElemIteratorPtr elemIt = smDS->GetElements();
546 if ( !elemIt->more() ) continue;
547 const SMDS_MeshElement* e = elemIt->next();
548 if ( !e || e->IsQuadratic() )
551 TIDSortedElemSet elems;
553 while ( elemIt->more() )
554 elems.insert( elems.end(), elemIt->next() );
556 SMESH_MeshEditor( mesh ).ConvertToQuadratic( /*3d=*/false, elems, /*biQuad=*/false );
562 //================================================================================
564 * \brief Initialize netgen::OCCGeometry with OCCT shape
566 //================================================================================
568 void NETGENPlugin_Mesher::PrepareOCCgeometry(netgen::OCCGeometry& occgeo,
569 const TopoDS_Shape& shape,
571 list< SMESH_subMesh* > * meshedSM,
572 NETGENPlugin_Internals* intern)
574 updateTriangulation( shape );
577 BRepBndLib::Add (shape, bb);
578 double x1,y1,z1,x2,y2,z2;
579 bb.Get (x1,y1,z1,x2,y2,z2);
580 MESSAGE("shape bounding box:\n" <<
581 "(" << x1 << " " << y1 << " " << z1 << ") " <<
582 "(" << x2 << " " << y2 << " " << z2 << ")");
583 netgen::Point<3> p1 = netgen::Point<3> (x1,y1,z1);
584 netgen::Point<3> p2 = netgen::Point<3> (x2,y2,z2);
585 occgeo.boundingbox = netgen::Box<3> (p1,p2);
587 occgeo.shape = shape;
590 // fill maps of shapes of occgeo with not yet meshed subshapes
592 // get root submeshes
593 list< SMESH_subMesh* > rootSM;
594 const int shapeID = mesh.GetMeshDS()->ShapeToIndex( shape );
595 if ( shapeID > 0 ) { // SMESH_subMesh with ID 0 may exist, don't use it!
596 rootSM.push_back( mesh.GetSubMesh( shape ));
599 for ( TopoDS_Iterator it( shape ); it.More(); it.Next() )
600 rootSM.push_back( mesh.GetSubMesh( it.Value() ));
605 // add subshapes of empty submeshes
606 list< SMESH_subMesh* >::iterator rootIt = rootSM.begin(), rootEnd = rootSM.end();
607 for ( ; rootIt != rootEnd; ++rootIt ) {
608 SMESH_subMesh * root = *rootIt;
609 SMESH_subMeshIteratorPtr smIt = root->getDependsOnIterator(/*includeSelf=*/true,
610 /*complexShapeFirst=*/true);
611 // to find a right orientation of subshapes (PAL20462)
612 TopTools_IndexedMapOfShape subShapes;
613 TopExp::MapShapes(root->GetSubShape(), subShapes);
614 while ( smIt->more() )
616 SMESH_subMesh* sm = smIt->next();
617 TopoDS_Shape shape = sm->GetSubShape();
618 totNbFaces += ( shape.ShapeType() == TopAbs_FACE );
619 if ( intern && intern->isShapeToPrecompute( shape ))
621 if ( !meshedSM || sm->IsEmpty() )
623 if ( shape.ShapeType() != TopAbs_VERTEX )
624 shape = subShapes( subShapes.FindIndex( shape ));// shape -> index -> oriented shape
625 if ( shape.Orientation() >= TopAbs_INTERNAL )
626 shape.Orientation( TopAbs_FORWARD ); // isuue 0020676
627 switch ( shape.ShapeType() ) {
628 case TopAbs_FACE : occgeo.fmap.Add( shape ); break;
629 case TopAbs_EDGE : occgeo.emap.Add( shape ); break;
630 case TopAbs_VERTEX: occgeo.vmap.Add( shape ); break;
631 case TopAbs_SOLID :occgeo.somap.Add( shape ); break;
635 // collect submeshes of meshed shapes
638 const int dim = SMESH_Gen::GetShapeDim( shape );
639 meshedSM[ dim ].push_back( sm );
643 occgeo.facemeshstatus.SetSize (totNbFaces);
644 occgeo.facemeshstatus = 0;
645 occgeo.face_maxh_modified.SetSize(totNbFaces);
646 occgeo.face_maxh_modified = 0;
647 occgeo.face_maxh.SetSize(totNbFaces);
648 occgeo.face_maxh = netgen::mparam.maxh;
651 //================================================================================
653 * \brief Return a default min size value suitable for the given geometry.
655 //================================================================================
657 double NETGENPlugin_Mesher::GetDefaultMinSize(const TopoDS_Shape& geom,
658 const double maxSize)
660 updateTriangulation( geom );
664 const int* pi[4] = { &i1, &i2, &i3, &i1 };
667 TopExp_Explorer fExp( geom, TopAbs_FACE );
668 for ( ; fExp.More(); fExp.Next() )
670 Handle(Poly_Triangulation) triangulation =
671 BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
672 if ( triangulation.IsNull() ) continue;
673 const double fTol = BRep_Tool::Tolerance( TopoDS::Face( fExp.Current() ));
674 const TColgp_Array1OfPnt& points = triangulation->Nodes();
675 const Poly_Array1OfTriangle& trias = triangulation->Triangles();
676 for ( int iT = trias.Lower(); iT <= trias.Upper(); ++iT )
678 trias(iT).Get( i1, i2, i3 );
679 for ( int j = 0; j < 3; ++j )
681 double dist2 = points(*pi[j]).SquareDistance( points( *pi[j+1] ));
682 if ( dist2 < minh && fTol*fTol < dist2 )
684 bb.Add( points(*pi[j]));
688 if ( minh > 0.25 * bb.SquareExtent() ) // simple geometry, rough triangulation
690 minh = 1e-3 * sqrt( bb.SquareExtent());
691 //cout << "BND BOX minh = " <<minh << endl;
695 minh = 3 * sqrt( minh ); // triangulation for visualization is rather fine
696 //cout << "TRIANGULATION minh = " <<minh << endl;
698 if ( minh > 0.5 * maxSize )
704 //================================================================================
706 * \brief Restrict size of elements at a given point
708 //================================================================================
710 void NETGENPlugin_Mesher::RestrictLocalSize(netgen::Mesh& ngMesh,
713 const bool overrideMinH)
715 if ( size <= std::numeric_limits<double>::min() )
717 if ( netgen::mparam.minh > size )
721 ngMesh.SetMinimalH( size );
722 netgen::mparam.minh = size;
726 size = netgen::mparam.minh;
729 netgen::Point3d pi(p.X(), p.Y(), p.Z());
730 ngMesh.RestrictLocalH( pi, size );
733 //================================================================================
735 * \brief fill ngMesh with nodes and elements of computed submeshes
737 //================================================================================
739 bool NETGENPlugin_Mesher::FillNgMesh(netgen::OCCGeometry& occgeom,
740 netgen::Mesh& ngMesh,
741 vector<const SMDS_MeshNode*>& nodeVec,
742 const list< SMESH_subMesh* > & meshedSM,
743 SMESH_MesherHelper* quadHelper,
744 SMESH_ProxyMesh::Ptr proxyMesh)
746 TNode2IdMap nodeNgIdMap;
747 for ( size_t i = 1; i < nodeVec.size(); ++i )
748 nodeNgIdMap.insert( make_pair( nodeVec[i], i ));
750 TopTools_MapOfShape visitedShapes;
751 map< SMESH_subMesh*, set< int > > visitedEdgeSM2Faces;
752 set< SMESH_subMesh* > computedSM( meshedSM.begin(), meshedSM.end() );
754 SMESH_MesherHelper helper (*_mesh);
756 int faceNgID = ngMesh.GetNFD();
758 list< SMESH_subMesh* >::const_iterator smIt, smEnd = meshedSM.end();
759 for ( smIt = meshedSM.begin(); smIt != smEnd; ++smIt )
761 SMESH_subMesh* sm = *smIt;
762 if ( !visitedShapes.Add( sm->GetSubShape() ))
765 const SMESHDS_SubMesh * smDS = sm->GetSubMeshDS();
766 if ( !smDS ) continue;
768 switch ( sm->GetSubShape().ShapeType() )
770 case TopAbs_EDGE: { // EDGE
771 // ----------------------
772 TopoDS_Edge geomEdge = TopoDS::Edge( sm->GetSubShape() );
773 if ( geomEdge.Orientation() >= TopAbs_INTERNAL )
774 geomEdge.Orientation( TopAbs_FORWARD ); // issue 0020676
776 // Add ng segments for each not meshed FACE the EDGE bounds
777 PShapeIteratorPtr fIt = helper.GetAncestors( geomEdge, *sm->GetFather(), TopAbs_FACE );
778 while ( const TopoDS_Shape * anc = fIt->next() )
780 faceNgID = occgeom.fmap.FindIndex( *anc );
782 continue; // meshed face
784 int faceSMDSId = helper.GetMeshDS()->ShapeToIndex( *anc );
785 if ( visitedEdgeSM2Faces[ sm ].count( faceSMDSId ))
786 continue; // already treated EDGE
788 TopoDS_Face face = TopoDS::Face( occgeom.fmap( faceNgID ));
789 if ( face.Orientation() >= TopAbs_INTERNAL )
790 face.Orientation( TopAbs_FORWARD ); // issue 0020676
792 // get all meshed EDGEs of the FACE connected to geomEdge (issue 0021140)
793 helper.SetSubShape( face );
794 list< TopoDS_Edge > edges = getConnectedEdges( geomEdge, face, computedSM, helper,
795 visitedEdgeSM2Faces );
797 continue; // wrong ancestor?
799 // find out orientation of <edges> within <face>
800 TopoDS_Edge eNotSeam = edges.front();
801 if ( helper.HasSeam() )
803 list< TopoDS_Edge >::iterator eIt = edges.begin();
804 while ( helper.IsRealSeam( *eIt )) ++eIt;
805 if ( eIt != edges.end() )
808 TopAbs_Orientation fOri = helper.GetSubShapeOri( face, eNotSeam );
809 bool isForwad = ( fOri == eNotSeam.Orientation() || fOri >= TopAbs_INTERNAL );
811 // get all nodes from connected <edges>
812 const bool isQuad = smDS->IsQuadratic();
813 StdMeshers_FaceSide fSide( face, edges, _mesh, isForwad, isQuad );
814 const vector<UVPtStruct>& points = fSide.GetUVPtStruct();
815 if ( points.empty() )
816 return false; // invalid node params?
817 int i, nbSeg = fSide.NbSegments();
819 // remember EDGEs of fSide to treat only once
820 for ( int iE = 0; iE < fSide.NbEdges(); ++iE )
821 visitedEdgeSM2Faces[ helper.GetMesh()->GetSubMesh( fSide.Edge(iE )) ].insert(faceSMDSId);
823 double otherSeamParam = 0;
828 int prevNgId = ngNodeId( points[0].node, ngMesh, nodeNgIdMap );
830 for ( i = 0; i < nbSeg; ++i )
832 const UVPtStruct& p1 = points[ i ];
833 const UVPtStruct& p2 = points[ i+1 ];
835 if ( p1.node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX ) //an EDGE begins
838 if ( helper.IsRealSeam( p1.node->getshapeId() ))
840 TopoDS_Edge e = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
841 isSeam = helper.IsRealSeam( e );
844 otherSeamParam = helper.GetOtherParam( helper.GetPeriodicIndex() & 1 ? p2.u : p2.v );
851 seg[1] = prevNgId = ngNodeId( p2.node, ngMesh, nodeNgIdMap );
852 // node param on curve
853 seg.epgeominfo[ 0 ].dist = p1.param;
854 seg.epgeominfo[ 1 ].dist = p2.param;
856 seg.epgeominfo[ 0 ].u = p1.u;
857 seg.epgeominfo[ 0 ].v = p1.v;
858 seg.epgeominfo[ 1 ].u = p2.u;
859 seg.epgeominfo[ 1 ].v = p2.v;
861 //geomEdge = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
862 //seg.epgeominfo[ 0 ].edgenr = seg.epgeominfo[ 1 ].edgenr = occgeom.emap.FindIndex( geomEdge );
864 //seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
865 seg.si = faceNgID; // = geom.fmap.FindIndex (face);
866 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
867 ngMesh.AddSegment (seg);
869 SMESH_TNodeXYZ np1( p1.node ), np2( p2.node );
870 RestrictLocalSize( ngMesh, 0.5*(np1+np2), (np1-np2).Modulus() );
873 cout << "Segment: " << seg.edgenr << " on SMESH face " << helper.GetMeshDS()->ShapeToIndex( face ) << endl
874 << "\tface index: " << seg.si << endl
875 << "\tp1: " << seg[0] << endl
876 << "\tp2: " << seg[1] << endl
877 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
878 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
879 //<< "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
880 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
881 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl;
882 //<< "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
886 if ( helper.GetPeriodicIndex() && 1 ) {
887 seg.epgeominfo[ 0 ].u = otherSeamParam;
888 seg.epgeominfo[ 1 ].u = otherSeamParam;
889 swap (seg.epgeominfo[0].v, seg.epgeominfo[1].v);
891 seg.epgeominfo[ 0 ].v = otherSeamParam;
892 seg.epgeominfo[ 1 ].v = otherSeamParam;
893 swap (seg.epgeominfo[0].u, seg.epgeominfo[1].u);
895 swap (seg[0], seg[1]);
896 swap (seg.epgeominfo[0].dist, seg.epgeominfo[1].dist);
897 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
898 ngMesh.AddSegment (seg);
900 cout << "Segment: " << seg.edgenr << endl
901 << "\t is SEAM (reverse) of the previous. "
902 << " Other " << (helper.GetPeriodicIndex() && 1 ? "U" : "V")
903 << " = " << otherSeamParam << endl;
906 else if ( fOri == TopAbs_INTERNAL )
908 swap (seg[0], seg[1]);
909 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
910 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
911 ngMesh.AddSegment (seg);
913 cout << "Segment: " << seg.edgenr << endl << "\t is REVERSE of the previous" << endl;
917 } // loop on geomEdge ancestors
919 if ( quadHelper ) // remember medium nodes of sub-meshes
921 SMDS_ElemIteratorPtr edges = smDS->GetElements();
922 while ( edges->more() )
924 const SMDS_MeshElement* e = edges->next();
925 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshEdge*>( e )))
931 } // case TopAbs_EDGE
933 case TopAbs_FACE: { // FACE
934 // ----------------------
935 const TopoDS_Face& geomFace = TopoDS::Face( sm->GetSubShape() );
936 helper.SetSubShape( geomFace );
937 bool isInternalFace = ( geomFace.Orientation() == TopAbs_INTERNAL );
939 // Find solids the geomFace bounds
940 int solidID1 = 0, solidID2 = 0;
941 StdMeshers_QuadToTriaAdaptor* quadAdaptor =
942 dynamic_cast<StdMeshers_QuadToTriaAdaptor*>( proxyMesh.get() );
945 solidID1 = occgeom.somap.FindIndex( quadAdaptor->GetShape() );
949 PShapeIteratorPtr solidIt = helper.GetAncestors( geomFace, *sm->GetFather(), TopAbs_SOLID);
950 while ( const TopoDS_Shape * solid = solidIt->next() )
952 int id = occgeom.somap.FindIndex ( *solid );
953 if ( solidID1 && id != solidID1 ) solidID2 = id;
957 // Add ng face descriptors of meshed faces
959 ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(faceNgID, solidID1, solidID2, 0));
961 // if second oreder is required, even already meshed faces must be passed to NETGEN
962 int fID = occgeom.fmap.Add( geomFace );
963 occgeom.facemeshstatus[ fID-1 ] = netgen::FACE_MESHED_OK;
964 while ( fID < faceNgID ) { // geomFace is already in occgeom.fmap, add a copy
965 fID = occgeom.fmap.Add( BRepBuilderAPI_Copy( geomFace, /*copyGeom=*/false ));
966 occgeom.facemeshstatus[ fID-1 ] = netgen::FACE_MESHED_OK;
968 // Problem with the second order in a quadrangular mesh remains.
969 // 1) All quadrangles generated by NETGEN are moved to an inexistent face
970 // by FillSMesh() (find "AddFaceDescriptor")
971 // 2) Temporary triangles generated by StdMeshers_QuadToTriaAdaptor
972 // are on faces where quadrangles were.
973 // Due to these 2 points, wrong geom faces are used while conversion to qudratic
974 // of the mentioned above quadrangles and triangles
976 // Orient the face correctly in solidID1 (issue 0020206)
977 bool reverse = false;
979 TopoDS_Shape solid = occgeom.somap( solidID1 );
980 TopAbs_Orientation faceOriInSolid = helper.GetSubShapeOri( solid, geomFace );
981 if ( faceOriInSolid >= 0 )
983 helper.IsReversedSubMesh( TopoDS::Face( geomFace.Oriented( faceOriInSolid )));
986 // Add surface elements
988 netgen::Element2d tri(3);
989 tri.SetIndex ( faceNgID );
990 SMESH_TNodeXYZ xyz[3];
992 #ifdef DUMP_TRIANGLES
993 cout << "SMESH face " << helper.GetMeshDS()->ShapeToIndex( geomFace )
994 << " internal="<<isInternalFace << endl;
997 smDS = proxyMesh->GetSubMesh( geomFace );
999 SMDS_ElemIteratorPtr faces = smDS->GetElements();
1000 while ( faces->more() )
1002 const SMDS_MeshElement* f = faces->next();
1003 if ( f->NbNodes() % 3 != 0 ) // not triangle
1005 PShapeIteratorPtr solidIt=helper.GetAncestors(geomFace,*sm->GetFather(),TopAbs_SOLID);
1006 if ( const TopoDS_Shape * solid = solidIt->next() )
1007 sm = _mesh->GetSubMesh( *solid );
1008 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1009 smError.reset( new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH,"Not triangle sub-mesh"));
1010 smError->myBadElements.push_back( f );
1014 for ( int i = 0; i < 3; ++i )
1016 const SMDS_MeshNode* node = f->GetNode( i ), * inFaceNode=0;
1019 // get node UV on face
1020 int shapeID = node->getshapeId();
1021 if ( helper.IsSeamShape( shapeID ))
1023 if ( helper.IsSeamShape( f->GetNodeWrap( i+1 )->getshapeId() ))
1024 inFaceNode = f->GetNodeWrap( i-1 );
1026 inFaceNode = f->GetNodeWrap( i+1 );
1028 gp_XY uv = helper.GetNodeUV( geomFace, node, inFaceNode );
1030 int ind = reverse ? 3-i : i+1;
1031 tri.GeomInfoPi(ind).u = uv.X();
1032 tri.GeomInfoPi(ind).v = uv.Y();
1033 tri.PNum (ind) = ngNodeId( node, ngMesh, nodeNgIdMap );
1036 // pass a triangle size to NG size-map
1037 double size = ( ( xyz[0] - xyz[1] ).Modulus() +
1038 ( xyz[1] - xyz[2] ).Modulus() +
1039 ( xyz[2] - xyz[0] ).Modulus() ) / 3;
1040 gp_XYZ gc = ( xyz[0] + xyz[1] + xyz[2] ) / 3;
1041 RestrictLocalSize( ngMesh, gc, size, /*overrideMinH=*/false );
1043 ngMesh.AddSurfaceElement (tri);
1044 #ifdef DUMP_TRIANGLES
1045 cout << tri << endl;
1048 if ( isInternalFace )
1050 swap( tri[1], tri[2] );
1051 ngMesh.AddSurfaceElement (tri);
1052 #ifdef DUMP_TRIANGLES
1053 cout << tri << endl;
1058 if ( quadHelper ) // remember medium nodes of sub-meshes
1060 SMDS_ElemIteratorPtr faces = smDS->GetElements();
1061 while ( faces->more() )
1063 const SMDS_MeshElement* f = faces->next();
1064 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshFace*>( f )))
1070 } // case TopAbs_FACE
1072 case TopAbs_VERTEX: { // VERTEX
1073 // --------------------------
1074 // issue 0021405. Add node only if a VERTEX is shared by a not meshed EDGE,
1075 // else netgen removes a free node and nodeVector becomes invalid
1076 PShapeIteratorPtr ansIt = helper.GetAncestors( sm->GetSubShape(),
1080 while ( const TopoDS_Shape* e = ansIt->next() )
1082 SMESH_subMesh* eSub = helper.GetMesh()->GetSubMesh( *e );
1083 if (( toAdd = ( eSub->IsEmpty() && !SMESH_Algo::isDegenerated( TopoDS::Edge( *e )))))
1088 SMDS_NodeIteratorPtr nodeIt = smDS->GetNodes();
1089 if ( nodeIt->more() )
1090 ngNodeId( nodeIt->next(), ngMesh, nodeNgIdMap );
1096 } // loop on submeshes
1099 nodeVec.resize( ngMesh.GetNP() + 1 );
1100 TNode2IdMap::iterator node_NgId, nodeNgIdEnd = nodeNgIdMap.end();
1101 for ( node_NgId = nodeNgIdMap.begin(); node_NgId != nodeNgIdEnd; ++node_NgId)
1102 nodeVec[ node_NgId->second ] = node_NgId->first;
1107 //================================================================================
1109 * \brief Duplicate mesh faces on internal geom faces
1111 //================================================================================
1113 void NETGENPlugin_Mesher::FixIntFaces(const netgen::OCCGeometry& occgeom,
1114 netgen::Mesh& ngMesh,
1115 NETGENPlugin_Internals& internalShapes)
1117 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1119 // find ng indices of internal faces
1121 for ( int ngFaceID = 1; ngFaceID <= occgeom.fmap.Extent(); ++ngFaceID )
1123 int smeshID = meshDS->ShapeToIndex( occgeom.fmap( ngFaceID ));
1124 if ( internalShapes.isInternalShape( smeshID ))
1125 ngFaceIds.insert( ngFaceID );
1127 if ( !ngFaceIds.empty() )
1130 int i, nbFaces = ngMesh.GetNSE();
1131 for ( i = 1; i <= nbFaces; ++i)
1133 netgen::Element2d elem = ngMesh.SurfaceElement(i);
1134 if ( ngFaceIds.count( elem.GetIndex() ))
1136 swap( elem[1], elem[2] );
1137 ngMesh.AddSurfaceElement (elem);
1143 //================================================================================
1145 * \brief Tries to heal the mesh on a FACE. The FACE is supposed to be partially
1146 * meshed due to NETGEN failure
1147 * \param [in] occgeom - geometry
1148 * \param [in,out] ngMesh - the mesh to fix
1149 * \param [inout] faceID - ID of the FACE to fix the mesh on
1150 * \return bool - is mesh is or becomes OK
1152 //================================================================================
1154 bool NETGENPlugin_Mesher::FixFaceMesh(const netgen::OCCGeometry& occgeom,
1155 netgen::Mesh& ngMesh,
1158 // we address a case where the FACE is almost fully meshed except small holes
1159 // of usually triangular shape at FACE boundary (IPAL52861)
1161 // The case appeared to be not simple: holes only look triangular but
1162 // indeed are a self intersecting polygon. A reason of the bug was in coincident
1163 // NG points on a seam edge. But the code below is very nice, leave it for
1168 if ( occgeom.fmap.Extent() < faceID )
1170 //const TopoDS_Face& face = TopoDS::Face( occgeom.fmap( faceID ));
1172 // find free links on the FACE
1173 NCollection_Map<Link> linkMap;
1174 for ( int iF = 1; iF <= ngMesh.GetNSE(); ++iF )
1176 const netgen::Element2d& elem = ngMesh.SurfaceElement(iF);
1177 if ( faceID != elem.GetIndex() )
1179 int n0 = elem[ elem.GetNP() - 1 ];
1180 for ( int i = 0; i < elem.GetNP(); ++i )
1183 Link link( n0, n1 );
1184 if ( !linkMap.Add( link ))
1185 linkMap.Remove( link );
1189 // add/remove boundary links
1190 for ( int iSeg = 1; iSeg <= ngMesh.GetNSeg(); ++iSeg )
1192 const netgen::Segment& seg = ngMesh.LineSegment( iSeg );
1193 if ( seg.si != faceID ) // !edgeIDs.Contains( seg.edgenr ))
1195 Link link( seg[1], seg[0] ); // reverse!!!
1196 if ( !linkMap.Add( link ))
1197 linkMap.Remove( link );
1199 if ( linkMap.IsEmpty() )
1201 if ( linkMap.Extent() < 3 )
1204 // make triangles of the links
1206 netgen::Element2d tri(3);
1207 tri.SetIndex ( faceID );
1209 NCollection_Map<Link>::Iterator linkIt( linkMap );
1210 Link link1 = linkIt.Value();
1211 // look for a link connected to link1
1212 NCollection_Map<Link>::Iterator linkIt2 = linkIt;
1213 for ( linkIt2.Next(); linkIt2.More(); linkIt2.Next() )
1215 const Link& link2 = linkIt2.Value();
1216 if ( link2.IsConnected( link1 ))
1218 // look for a link connected to both link1 and link2
1219 NCollection_Map<Link>::Iterator linkIt3 = linkIt2;
1220 for ( linkIt3.Next(); linkIt3.More(); linkIt3.Next() )
1222 const Link& link3 = linkIt3.Value();
1223 if ( link3.IsConnected( link1 ) &&
1224 link3.IsConnected( link2 ) )
1229 tri[2] = ( link2.Contains( link1.n1 ) ? link2.n1 : link3.n1 );
1230 if ( tri[0] == tri[2] || tri[1] == tri[2] )
1232 ngMesh.AddSurfaceElement( tri );
1234 // prepare for the next tria search
1235 if ( linkMap.Extent() == 3 )
1237 linkMap.Remove( link3 );
1238 linkMap.Remove( link2 );
1240 linkMap.Remove( link1 );
1241 link1 = linkIt.Value();
1254 //================================================================================
1255 // define gp_XY_Subtracted pointer to function calling gp_XY::Subtracted(gp_XY)
1256 gp_XY_FunPtr(Subtracted);
1257 //gp_XY_FunPtr(Added);
1259 //================================================================================
1261 * \brief Evaluate distance between two 2d points along the surface
1263 //================================================================================
1265 double evalDist( const gp_XY& uv1,
1267 const Handle(Geom_Surface)& surf,
1268 const int stopHandler=-1)
1270 if ( stopHandler > 0 ) // continue recursion
1272 gp_XY mid = SMESH_MesherHelper::GetMiddleUV( surf, uv1, uv2 );
1273 return evalDist( uv1,mid, surf, stopHandler-1 ) + evalDist( mid,uv2, surf, stopHandler-1 );
1275 double dist3D = surf->Value( uv1.X(), uv1.Y() ).Distance( surf->Value( uv2.X(), uv2.Y() ));
1276 if ( stopHandler == 0 ) // stop recursion
1279 // start recursion if necessary
1280 double dist2D = SMESH_MesherHelper::ApplyIn2D(surf, uv1, uv2, gp_XY_Subtracted, 0).Modulus();
1281 if ( fabs( dist3D - dist2D ) < dist2D * 1e-10 )
1282 return dist3D; // equal parametrization of a planar surface
1284 return evalDist( uv1, uv2, surf, 3 ); // start recursion
1287 //================================================================================
1289 * \brief Data of vertex internal in geom face
1291 //================================================================================
1295 gp_XY uv; //!< UV in face parametric space
1296 int ngId; //!< ng id of corrsponding node
1297 gp_XY uvClose; //!< UV of closest boundary node
1298 int ngIdClose; //!< ng id of closest boundary node
1301 //================================================================================
1303 * \brief Data of vertex internal in solid
1305 //================================================================================
1309 int ngId; //!< ng id of corresponding node
1310 int ngIdClose; //!< ng id of closest 2d mesh element
1311 int ngIdCloseN; //!< ng id of closest node of the closest 2d mesh element
1314 inline double dist2(const netgen::MeshPoint& p1, const netgen::MeshPoint& p2)
1316 return gp_Pnt( NGPOINT_COORDS(p1)).SquareDistance( gp_Pnt( NGPOINT_COORDS(p2)));
1320 //================================================================================
1322 * \brief Make netgen take internal vertices in faces into account by adding
1323 * segments including internal vertices
1325 * This function works in supposition that 1D mesh is already computed in ngMesh
1327 //================================================================================
1329 void NETGENPlugin_Mesher::AddIntVerticesInFaces(const netgen::OCCGeometry& occgeom,
1330 netgen::Mesh& ngMesh,
1331 vector<const SMDS_MeshNode*>& nodeVec,
1332 NETGENPlugin_Internals& internalShapes)
1334 if ((int) nodeVec.size() < ngMesh.GetNP() )
1335 nodeVec.resize( ngMesh.GetNP(), 0 );
1337 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1338 SMESH_MesherHelper helper( internalShapes.getMesh() );
1340 const map<int,list<int> >& face2Vert = internalShapes.getFacesWithVertices();
1341 map<int,list<int> >::const_iterator f2v = face2Vert.begin();
1342 for ( ; f2v != face2Vert.end(); ++f2v )
1344 const TopoDS_Face& face = TopoDS::Face( meshDS->IndexToShape( f2v->first ));
1345 if ( face.IsNull() ) continue;
1346 int faceNgID = occgeom.fmap.FindIndex (face);
1347 if ( faceNgID < 0 ) continue;
1349 TopLoc_Location loc;
1350 Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc);
1352 helper.SetSubShape( face );
1353 helper.SetElementsOnShape( true );
1355 // Get data of internal vertices and add them to ngMesh
1357 multimap< double, TIntVData > dist2VData; // sort vertices by distance from boundary nodes
1359 int i, nbSegInit = ngMesh.GetNSeg();
1361 // boundary characteristics
1362 double totSegLen2D = 0;
1365 const list<int>& iVertices = f2v->second;
1366 list<int>::const_iterator iv = iVertices.begin();
1367 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1370 // get node on vertex
1371 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1372 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1375 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1376 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1377 nV = SMESH_Algo::VertexNode( V, meshDS );
1378 if ( !nV ) continue;
1381 netgen::MeshPoint mp( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1382 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1383 vData.ngId = ngMesh.GetNP();
1384 nodeVec.push_back( nV );
1388 vData.uv = helper.GetNodeUV( face, nV, 0, &uvOK );
1389 if ( !uvOK ) helper.CheckNodeUV( face, nV, vData.uv, BRep_Tool::Tolerance(V),/*force=*/1);
1391 // loop on all segments of the face to find the node closest to vertex and to count
1392 // average segment 2d length
1393 double closeDist2 = numeric_limits<double>::max(), dist2;
1395 for (i = 1; i <= ngMesh.GetNSeg(); ++i)
1397 netgen::Segment & seg = ngMesh.LineSegment(i);
1398 if ( seg.si != faceNgID ) continue;
1400 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1402 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1403 if ( ngIdLast == seg[ iEnd ] ) continue;
1404 dist2 = helper.ApplyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1405 if ( dist2 < closeDist2 )
1406 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1407 ngIdLast = seg[ iEnd ];
1411 totSegLen2D += helper.ApplyIn2D(surf, uv[0], uv[1], gp_XY_Subtracted, false).Modulus();
1415 dist2VData.insert( make_pair( closeDist2, vData ));
1418 if ( totNbSeg == 0 ) break;
1419 double avgSegLen2d = totSegLen2D / totNbSeg;
1421 // Loop on vertices to add segments
1423 multimap< double, TIntVData >::iterator dist_vData = dist2VData.begin();
1424 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1426 double closeDist2 = dist_vData->first, dist2;
1427 TIntVData & vData = dist_vData->second;
1429 // try to find more close node among segments added for internal vertices
1430 for (i = nbSegInit+1; i <= ngMesh.GetNSeg(); ++i)
1432 netgen::Segment & seg = ngMesh.LineSegment(i);
1433 if ( seg.si != faceNgID ) continue;
1435 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1437 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1438 dist2 = helper.ApplyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1439 if ( dist2 < closeDist2 )
1440 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1443 // decide whether to use the closest node as the second end of segment or to
1444 // create a new point
1445 int segEnd1 = vData.ngId;
1446 int segEnd2 = vData.ngIdClose; // to use closest node
1447 gp_XY uvV = vData.uv, uvP = vData.uvClose;
1448 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1449 double nodeDist2D = sqrt( closeDist2 );
1450 double nodeDist3D = evalDist( vData.uv, vData.uvClose, surf );
1451 bool avgLenOK = ( avgSegLen2d < 0.75 * nodeDist2D );
1452 bool hintLenOK = ( segLenHint < 0.75 * nodeDist3D );
1453 //cout << "uvV " << uvV.X() <<","<<uvV.Y() << " ";
1454 if ( hintLenOK || avgLenOK )
1456 // create a point between the closest node and V
1459 double r = min( 0.5, ( hintLenOK ? segLenHint/nodeDist3D : avgSegLen2d/nodeDist2D ));
1460 // direction from V to closet node in 2D
1461 gp_Dir2d v2n( helper.ApplyIn2D(surf, uvP, uvV, gp_XY_Subtracted, false ));
1463 uvP = vData.uv + r * nodeDist2D * v2n.XY();
1464 gp_Pnt P = surf->Value( uvP.X(), uvP.Y() ).Transformed( loc );
1466 netgen::MeshPoint mp( netgen::Point<3> (P.X(), P.Y(), P.Z()));
1467 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1468 segEnd2 = ngMesh.GetNP();
1469 //cout << "Middle " << r << " uv " << uvP.X() << "," << uvP.Y() << "( " << ngMesh.Point(segEnd2).X()<<","<<ngMesh.Point(segEnd2).Y()<<","<<ngMesh.Point(segEnd2).Z()<<" )"<< endl;
1470 SMDS_MeshNode * nP = helper.AddNode(P.X(), P.Y(), P.Z());
1471 nodeVec.push_back( nP );
1473 //else cout << "at Node " << " uv " << uvP.X() << "," << uvP.Y() << endl;
1476 netgen::Segment seg;
1478 if ( segEnd1 > segEnd2 ) swap( segEnd1, segEnd2 ), swap( uvV, uvP );
1479 seg[0] = segEnd1; // ng node id
1480 seg[1] = segEnd2; // ng node id
1481 seg.edgenr = ngMesh.GetNSeg() + 1;// segment id
1484 seg.epgeominfo[ 0 ].dist = 0; // param on curve
1485 seg.epgeominfo[ 0 ].u = uvV.X();
1486 seg.epgeominfo[ 0 ].v = uvV.Y();
1487 seg.epgeominfo[ 1 ].dist = 1; // param on curve
1488 seg.epgeominfo[ 1 ].u = uvP.X();
1489 seg.epgeominfo[ 1 ].v = uvP.Y();
1491 // seg.epgeominfo[ 0 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1492 // seg.epgeominfo[ 1 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1494 ngMesh.AddSegment (seg);
1496 // add reverse segment
1497 swap (seg[0], seg[1]);
1498 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1499 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1500 ngMesh.AddSegment (seg);
1506 //================================================================================
1508 * \brief Make netgen take internal vertices in solids into account by adding
1509 * faces including internal vertices
1511 * This function works in supposition that 2D mesh is already computed in ngMesh
1513 //================================================================================
1515 void NETGENPlugin_Mesher::AddIntVerticesInSolids(const netgen::OCCGeometry& occgeom,
1516 netgen::Mesh& ngMesh,
1517 vector<const SMDS_MeshNode*>& nodeVec,
1518 NETGENPlugin_Internals& internalShapes)
1520 #ifdef DUMP_TRIANGLES_SCRIPT
1521 // create a python script making a mesh containing triangles added for internal vertices
1522 ofstream py(DUMP_TRIANGLES_SCRIPT);
1523 py << "import SMESH"<< endl
1524 << "from salome.smesh import smeshBuilder"<<endl
1525 << "smesh = smeshBuilder.New(salome.myStudy)"<<endl
1526 << "m = smesh.Mesh(name='triangles')" << endl;
1528 if ((int) nodeVec.size() < ngMesh.GetNP() )
1529 nodeVec.resize( ngMesh.GetNP(), 0 );
1531 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1532 SMESH_MesherHelper helper( internalShapes.getMesh() );
1534 const map<int,list<int> >& so2Vert = internalShapes.getSolidsWithVertices();
1535 map<int,list<int> >::const_iterator s2v = so2Vert.begin();
1536 for ( ; s2v != so2Vert.end(); ++s2v )
1538 const TopoDS_Shape& solid = meshDS->IndexToShape( s2v->first );
1539 if ( solid.IsNull() ) continue;
1540 int solidNgID = occgeom.somap.FindIndex (solid);
1541 if ( solidNgID < 0 && !occgeom.somap.IsEmpty() ) continue;
1543 helper.SetSubShape( solid );
1544 helper.SetElementsOnShape( true );
1546 // find ng indices of faces within the solid
1548 for (TopExp_Explorer fExp(solid, TopAbs_FACE); fExp.More(); fExp.Next() )
1549 ngFaceIds.insert( occgeom.fmap.FindIndex( fExp.Current() ));
1550 if ( ngFaceIds.size() == 1 && *ngFaceIds.begin() == 0 )
1551 ngFaceIds.insert( 1 );
1553 // Get data of internal vertices and add them to ngMesh
1555 multimap< double, TIntVSoData > dist2VData; // sort vertices by distance from ng faces
1557 int i, nbFaceInit = ngMesh.GetNSE();
1559 // boundary characteristics
1560 double totSegLen = 0;
1563 const list<int>& iVertices = s2v->second;
1564 list<int>::const_iterator iv = iVertices.begin();
1565 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1568 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1570 // get node on vertex
1571 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1574 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1575 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1576 nV = SMESH_Algo::VertexNode( V, meshDS );
1577 if ( !nV ) continue;
1580 netgen::MeshPoint mpV( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1581 ngMesh.AddPoint ( mpV, 1, netgen::FIXEDPOINT );
1582 vData.ngId = ngMesh.GetNP();
1583 nodeVec.push_back( nV );
1585 // loop on all 2d elements to find the one closest to vertex and to count
1586 // average segment length
1587 double closeDist2 = numeric_limits<double>::max(), avgDist2;
1588 for (i = 1; i <= ngMesh.GetNSE(); ++i)
1590 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1591 if ( !ngFaceIds.count( elem.GetIndex() )) continue;
1593 multimap< double, int> dist2nID; // sort nodes of element by distance from V
1594 for ( int j = 0; j < elem.GetNP(); ++j)
1596 netgen::MeshPoint mp = ngMesh.Point( elem[j] );
1597 double d2 = dist2( mpV, mp );
1598 dist2nID.insert( make_pair( d2, elem[j] ));
1599 avgDist2 += d2 / elem.GetNP();
1601 totNbSeg++, totSegLen+= sqrt( dist2( mp, ngMesh.Point( elem[(j+1)%elem.GetNP()])));
1603 double dist = dist2nID.begin()->first; //avgDist2;
1604 if ( dist < closeDist2 )
1605 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= dist;
1607 dist2VData.insert( make_pair( closeDist2, vData ));
1610 if ( totNbSeg == 0 ) break;
1611 double avgSegLen = totSegLen / totNbSeg;
1613 // Loop on vertices to add triangles
1615 multimap< double, TIntVSoData >::iterator dist_vData = dist2VData.begin();
1616 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1618 double closeDist2 = dist_vData->first;
1619 TIntVSoData & vData = dist_vData->second;
1621 const netgen::MeshPoint& mpV = ngMesh.Point( vData.ngId );
1623 // try to find more close face among ones added for internal vertices
1624 for (i = nbFaceInit+1; i <= ngMesh.GetNSE(); ++i)
1626 double avgDist2 = 0;
1627 multimap< double, int> dist2nID;
1628 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1629 for ( int j = 0; j < elem.GetNP(); ++j)
1631 double d = dist2( mpV, ngMesh.Point( elem[j] ));
1632 dist2nID.insert( make_pair( d, elem[j] ));
1633 avgDist2 += d / elem.GetNP();
1634 if ( avgDist2 < closeDist2 )
1635 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= avgDist2;
1638 // sort nodes of the closest face by angle with vector from V to the closest node
1639 const double tol = numeric_limits<double>::min();
1640 map< double, int > angle2ID;
1641 const netgen::Element2d& closeFace = ngMesh.SurfaceElement( vData.ngIdClose );
1642 netgen::MeshPoint mp[2];
1643 mp[0] = ngMesh.Point( vData.ngIdCloseN );
1644 gp_XYZ p1( NGPOINT_COORDS( mp[0] ));
1645 gp_XYZ pV( NGPOINT_COORDS( mpV ));
1646 gp_Vec v2p1( pV, p1 );
1647 double distN1 = v2p1.Magnitude();
1648 if ( distN1 <= tol ) continue;
1650 for ( int j = 0; j < closeFace.GetNP(); ++j)
1652 mp[1] = ngMesh.Point( closeFace[j] );
1653 gp_Vec v2p( pV, gp_Pnt( NGPOINT_COORDS( mp[1] )) );
1654 angle2ID.insert( make_pair( v2p1.Angle( v2p ), closeFace[j]));
1656 // get node with angle of 60 degrees or greater
1657 map< double, int >::iterator angle_id = angle2ID.lower_bound( 60. * M_PI / 180. );
1658 if ( angle_id == angle2ID.end() ) angle_id = --angle2ID.end();
1659 const double minAngle = 30. * M_PI / 180.;
1660 const double angle = angle_id->first;
1661 bool angleOK = ( angle > minAngle );
1663 // find points to create a triangle
1664 netgen::Element2d tri(3);
1666 tri[0] = vData.ngId;
1667 tri[1] = vData.ngIdCloseN; // to use the closest nodes
1668 tri[2] = angle_id->second; // to use the node with best angle
1670 // decide whether to use the closest node and the node with best angle or to create new ones
1671 for ( int isBestAngleN = 0; isBestAngleN < 2; ++isBestAngleN )
1673 bool createNew = !angleOK; //, distOK = true;
1675 int triInd = isBestAngleN ? 2 : 1;
1676 mp[isBestAngleN] = ngMesh.Point( tri[triInd] );
1681 double distN2 = sqrt( dist2( mpV, mp[isBestAngleN]));
1682 createNew = ( fabs( distN2 - distN1 ) > 0.25 * distN1 );
1684 else if ( angle < tol )
1686 v2p1.SetX( v2p1.X() + 1e-3 );
1692 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1693 bool avgLenOK = ( avgSegLen < 0.75 * distN1 );
1694 bool hintLenOK = ( segLenHint < 0.75 * distN1 );
1695 createNew = (createNew || avgLenOK || hintLenOK );
1696 // we create a new node not closer than 0.5 to the closest face
1697 // in order not to clash with other close face
1698 double r = min( 0.5, ( hintLenOK ? segLenHint : avgSegLen ) / distN1 );
1699 distFromV = r * distN1;
1703 // create a new point, between the node and the vertex if angleOK
1704 gp_XYZ p( NGPOINT_COORDS( mp[isBestAngleN] ));
1705 gp_Vec v2p( pV, p ); v2p.Normalize();
1706 if ( isBestAngleN && !angleOK )
1707 p = p1 + gp_Dir( v2p.XYZ() - v2p1.XYZ()).XYZ() * distN1 * 0.95;
1709 p = pV + v2p.XYZ() * distFromV;
1711 if ( !isBestAngleN ) p1 = p, distN1 = distFromV;
1713 mp[isBestAngleN].SetPoint( netgen::Point<3> (p.X(), p.Y(), p.Z()));
1714 ngMesh.AddPoint ( mp[isBestAngleN], 1, netgen::SURFACEPOINT );
1715 tri[triInd] = ngMesh.GetNP();
1716 nodeVec.push_back( helper.AddNode( p.X(), p.Y(), p.Z()) );
1719 ngMesh.AddSurfaceElement (tri);
1720 swap( tri[1], tri[2] );
1721 ngMesh.AddSurfaceElement (tri);
1723 #ifdef DUMP_TRIANGLES_SCRIPT
1724 py << "n1 = m.AddNode( "<< mpV(0)<<", "<< mpV(1)<<", "<< mpV(2)<<") "<< endl
1725 << "n2 = m.AddNode( "<< mp[0](0)<<", "<< mp[0](1)<<", "<< mp[0](2)<<") "<< endl
1726 << "n3 = m.AddNode( "<< mp[1](0)<<", "<< mp[1](1)<<", "<< mp[1](2)<<" )" << endl
1727 << "m.AddFace([n1,n2,n3])" << endl;
1729 } // loop on internal vertices of a solid
1731 } // loop on solids with internal vertices
1734 //================================================================================
1736 * \brief Fill netgen mesh with segments of a FACE
1737 * \param ngMesh - netgen mesh
1738 * \param geom - container of OCCT geometry to mesh
1739 * \param wires - data of nodes on FACE boundary
1740 * \param helper - mesher helper holding the FACE
1741 * \param nodeVec - vector of nodes in which node index == netgen ID
1742 * \retval SMESH_ComputeErrorPtr - error description
1744 //================================================================================
1746 SMESH_ComputeErrorPtr
1747 NETGENPlugin_Mesher::AddSegmentsToMesh(netgen::Mesh& ngMesh,
1748 netgen::OCCGeometry& geom,
1749 const TSideVector& wires,
1750 SMESH_MesherHelper& helper,
1751 vector< const SMDS_MeshNode* > & nodeVec,
1752 const bool overrideMinH)
1754 // ----------------------------
1755 // Check wires and count nodes
1756 // ----------------------------
1758 for ( size_t iW = 0; iW < wires.size(); ++iW )
1760 StdMeshers_FaceSidePtr wire = wires[ iW ];
1761 if ( wire->MissVertexNode() )
1763 // Commented for issue 0020960. It worked for the case, let's wait for case where it doesn't.
1764 // It seems that there is no reason for this limitation
1766 // (new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH, "Missing nodes on vertices"));
1768 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1769 if ((int) uvPtVec.size() != wire->NbPoints() )
1770 return SMESH_ComputeError::New(COMPERR_BAD_INPUT_MESH,
1771 SMESH_Comment("Unexpected nb of points on wire ") << iW
1772 << ": " << uvPtVec.size()<<" != "<<wire->NbPoints());
1773 nbNodes += wire->NbPoints();
1775 nodeVec.reserve( nodeVec.size() + nbNodes + 1 );
1776 if ( nodeVec.empty() )
1777 nodeVec.push_back( 0 );
1779 // -----------------
1781 // -----------------
1783 const bool wasNgMeshEmpty = ( ngMesh.GetNP() < 1 ); /* true => this method is called by
1784 NETGENPlugin_NETGEN_2D_ONLY */
1786 // map for nodes on vertices since they can be shared between wires
1787 // ( issue 0020676, face_int_box.brep) and nodes built by NETGEN
1788 map<const SMDS_MeshNode*, int > node2ngID;
1789 if ( !wasNgMeshEmpty ) // fill node2ngID with nodes built by NETGEN
1791 set< int > subIDs; // ids of sub-shapes of the FACE
1792 for ( size_t iW = 0; iW < wires.size(); ++iW )
1794 StdMeshers_FaceSidePtr wire = wires[ iW ];
1795 for ( int iE = 0, nbE = wire->NbEdges(); iE < nbE; ++iE )
1797 subIDs.insert( wire->EdgeID( iE ));
1798 subIDs.insert( helper.GetMeshDS()->ShapeToIndex( wire->FirstVertex( iE )));
1801 for ( size_t ngID = 1; ngID < nodeVec.size(); ++ngID )
1802 if ( subIDs.count( nodeVec[ngID]->getshapeId() ))
1803 node2ngID.insert( make_pair( nodeVec[ngID], ngID ));
1806 const int solidID = 0, faceID = geom.fmap.FindIndex( helper.GetSubShape() );
1807 if ( ngMesh.GetNFD() < 1 )
1808 ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(faceID, solidID, solidID, 0));
1810 for ( size_t iW = 0; iW < wires.size(); ++iW )
1812 StdMeshers_FaceSidePtr wire = wires[ iW ];
1813 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1814 const int nbSegments = wire->NbPoints() - 1;
1816 // assure the 1st node to be in node2ngID, which is needed to correctly
1817 // "close chain of segments" (see below) in case if the 1st node is not
1818 // onVertex because it is on a Viscous layer
1819 node2ngID.insert( make_pair( uvPtVec[ 0 ].node, ngMesh.GetNP() + 1 ));
1821 // compute length of every segment
1822 vector<double> segLen( nbSegments );
1823 for ( int i = 0; i < nbSegments; ++i )
1824 segLen[i] = SMESH_TNodeXYZ( uvPtVec[ i ].node ).Distance( uvPtVec[ i+1 ].node );
1826 int edgeID = 1, posID = -2;
1827 bool isInternalWire = false;
1828 double vertexNormPar = 0;
1829 //const int prevNbNGSeg = ngMesh.GetNSeg();
1830 for ( int i = 0; i < nbSegments; ++i ) // loop on segments
1832 // Add the first point of a segment
1834 const SMDS_MeshNode * n = uvPtVec[ i ].node;
1835 const int posShapeID = n->getshapeId();
1836 bool onVertex = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX );
1837 bool onEdge = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_EDGE );
1839 // skip nodes on degenerated edges
1840 if ( helper.IsDegenShape( posShapeID ) &&
1841 helper.IsDegenShape( uvPtVec[ i+1 ].node->getshapeId() ))
1844 int ngID1 = ngMesh.GetNP() + 1, ngID2 = ngID1+1;
1845 if ( onVertex || ( !wasNgMeshEmpty && onEdge ) || helper.IsRealSeam( posShapeID ))
1846 ngID1 = node2ngID.insert( make_pair( n, ngID1 )).first->second;
1847 if ( ngID1 > ngMesh.GetNP() )
1849 netgen::MeshPoint mp( netgen::Point<3> (n->X(), n->Y(), n->Z()) );
1850 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1851 nodeVec.push_back( n );
1853 else // n is in ngMesh already, and ngID2 in prev segment is wrong
1855 ngID2 = ngMesh.GetNP() + 1;
1856 if ( i > 0 ) // prev segment belongs to same wire
1858 netgen::Segment& prevSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
1865 netgen::Segment seg;
1867 seg[0] = ngID1; // ng node id
1868 seg[1] = ngID2; // ng node id
1869 seg.edgenr = ngMesh.GetNSeg() + 1; // ng segment id
1870 seg.si = faceID; // = geom.fmap.FindIndex (face);
1872 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1874 const UVPtStruct& pnt = uvPtVec[ i + iEnd ];
1876 seg.epgeominfo[ iEnd ].dist = pnt.param; // param on curve
1877 seg.epgeominfo[ iEnd ].u = pnt.u;
1878 seg.epgeominfo[ iEnd ].v = pnt.v;
1880 // find out edge id and node parameter on edge
1881 onVertex = ( pnt.normParam + 1e-10 > vertexNormPar );
1882 if ( onVertex || posShapeID != posID )
1885 double normParam = pnt.normParam;
1887 normParam = 0.5 * ( uvPtVec[ i ].normParam + uvPtVec[ i+1 ].normParam );
1888 int edgeIndexInWire = wire->EdgeIndex( normParam );
1889 vertexNormPar = wire->LastParameter( edgeIndexInWire );
1890 const TopoDS_Edge& edge = wire->Edge( edgeIndexInWire );
1891 edgeID = geom.emap.FindIndex( edge );
1893 isInternalWire = ( edge.Orientation() == TopAbs_INTERNAL );
1894 // if ( onVertex ) // param on curve is different on each of two edges
1895 // seg.epgeominfo[ iEnd ].dist = helper.GetNodeU( edge, pnt.node );
1897 seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
1900 ngMesh.AddSegment (seg);
1902 // restrict size of elements near the segment
1903 SMESH_TNodeXYZ np1( n ), np2( uvPtVec[ i+1 ].node );
1904 // get an average size of adjacent segments to avoid sharp change of
1905 // element size (regression on issue 0020452, note 0010898)
1906 int iPrev = SMESH_MesherHelper::WrapIndex( i-1, nbSegments );
1907 int iNext = SMESH_MesherHelper::WrapIndex( i+1, nbSegments );
1908 double sumH = segLen[ iPrev ] + segLen[ i ] + segLen[ iNext ];
1909 int nbSeg = ( int( segLen[ iPrev ] > sumH / 100.) +
1910 int( segLen[ i ] > sumH / 100.) +
1911 int( segLen[ iNext ] > sumH / 100.));
1913 RestrictLocalSize( ngMesh, 0.5*(np1+np2), sumH / nbSeg, overrideMinH );
1915 if ( isInternalWire )
1917 swap (seg[0], seg[1]);
1918 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1919 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1920 ngMesh.AddSegment (seg);
1922 } // loop on segments on a wire
1924 // close chain of segments
1925 if ( nbSegments > 0 )
1927 netgen::Segment& lastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() - int( isInternalWire));
1928 const SMDS_MeshNode * lastNode = uvPtVec.back().node;
1929 lastSeg[1] = node2ngID.insert( make_pair( lastNode, lastSeg[1] )).first->second;
1930 if ( lastSeg[1] > ngMesh.GetNP() )
1932 netgen::MeshPoint mp( netgen::Point<3> (lastNode->X(), lastNode->Y(), lastNode->Z()) );
1933 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1934 nodeVec.push_back( lastNode );
1936 if ( isInternalWire )
1938 netgen::Segment& realLastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
1939 realLastSeg[0] = lastSeg[1];
1943 #ifdef DUMP_SEGMENTS
1944 cout << "BEGIN WIRE " << iW << endl;
1945 for ( int i = prevNbNGSeg+1; i <= ngMesh.GetNSeg(); ++i )
1947 netgen::Segment& seg = ngMesh.LineSegment( i );
1949 netgen::Segment& prevSeg = ngMesh.LineSegment( i-1 );
1950 if ( seg[0] == prevSeg[1] && seg[1] == prevSeg[0] )
1952 cout << "Segment: " << seg.edgenr << endl << "\tis REVRESE of the previous one" << endl;
1956 cout << "Segment: " << seg.edgenr << endl
1957 << "\tp1: " << seg[0] << " n" << nodeVec[ seg[0]]->GetID() << endl
1958 << "\tp2: " << seg[1] << " n" << nodeVec[ seg[1]]->GetID() << endl
1959 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
1960 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
1961 << "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
1962 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
1963 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl
1964 << "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
1966 cout << "--END WIRE " << iW << endl;
1969 } // loop on WIREs of a FACE
1971 // add a segment instead of an internal vertex
1972 if ( wasNgMeshEmpty )
1974 NETGENPlugin_Internals intShapes( *helper.GetMesh(), helper.GetSubShape(), /*is3D=*/false );
1975 AddIntVerticesInFaces( geom, ngMesh, nodeVec, intShapes );
1977 ngMesh.CalcSurfacesOfNode();
1982 //================================================================================
1984 * \brief Fill SMESH mesh according to contents of netgen mesh
1985 * \param occgeo - container of OCCT geometry to mesh
1986 * \param ngMesh - netgen mesh
1987 * \param initState - bn of entities in netgen mesh before computing
1988 * \param sMesh - SMESH mesh to fill in
1989 * \param nodeVec - vector of nodes in which node index == netgen ID
1990 * \param comment - returns problem description
1991 * \param quadHelper - holder of medium nodes of sub-meshes
1992 * \retval int - error
1994 //================================================================================
1996 int NETGENPlugin_Mesher::FillSMesh(const netgen::OCCGeometry& occgeo,
1997 netgen::Mesh& ngMesh,
1998 const NETGENPlugin_ngMeshInfo& initState,
2000 std::vector<const SMDS_MeshNode*>& nodeVec,
2001 SMESH_Comment& comment,
2002 SMESH_MesherHelper* quadHelper)
2004 int nbNod = ngMesh.GetNP();
2005 int nbSeg = ngMesh.GetNSeg();
2006 int nbFac = ngMesh.GetNSE();
2007 int nbVol = ngMesh.GetNE();
2009 SMESHDS_Mesh* meshDS = sMesh.GetMeshDS();
2011 // quadHelper is used for either
2012 // 1) making quadratic elements when a lower dimention mesh is loaded
2013 // to SMESH before convertion to quadratic by NETGEN
2014 // 2) sewing of quadratic elements with quadratic elements of sub-meshes
2015 if ( quadHelper && !quadHelper->GetIsQuadratic() && quadHelper->GetTLinkNodeMap().empty() )
2018 // -------------------------------------
2019 // Create and insert nodes into nodeVec
2020 // -------------------------------------
2022 nodeVec.resize( nbNod + 1 );
2023 int i, nbInitNod = initState._nbNodes;
2024 for (i = nbInitNod+1; i <= nbNod; ++i )
2026 const netgen::MeshPoint& ngPoint = ngMesh.Point(i);
2027 SMDS_MeshNode* node = NULL;
2028 TopoDS_Vertex aVert;
2029 // First, netgen creates nodes on vertices in occgeo.vmap,
2030 // so node index corresponds to vertex index
2031 // but (issue 0020776) netgen does not create nodes with equal coordinates
2032 if ( i-nbInitNod <= occgeo.vmap.Extent() )
2034 gp_Pnt p ( NGPOINT_COORDS(ngPoint) );
2035 for (int iV = i-nbInitNod; aVert.IsNull() && iV <= occgeo.vmap.Extent(); ++iV)
2037 aVert = TopoDS::Vertex( occgeo.vmap( iV ));
2038 gp_Pnt pV = BRep_Tool::Pnt( aVert );
2039 if ( p.SquareDistance( pV ) > 1e-20 )
2042 node = const_cast<SMDS_MeshNode*>( SMESH_Algo::VertexNode( aVert, meshDS ));
2045 if (!node) // node not found on vertex
2047 node = meshDS->AddNode( NGPOINT_COORDS( ngPoint ));
2048 if (!aVert.IsNull())
2049 meshDS->SetNodeOnVertex(node, aVert);
2054 // -------------------------------------------
2055 // Create mesh segments along geometric edges
2056 // -------------------------------------------
2058 int nbInitSeg = initState._nbSegments;
2059 for (i = nbInitSeg+1; i <= nbSeg; ++i )
2061 const netgen::Segment& seg = ngMesh.LineSegment(i);
2063 int pinds[3] = { seg.pnums[0], seg.pnums[1], seg.pnums[2] };
2066 for (int j=0; j < 3; ++j)
2068 int pind = pinds[j];
2069 if (pind <= 0 || !nodeVec_ACCESS(pind))
2077 int aGeomEdgeInd = seg.epgeominfo[j].edgenr;
2078 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
2079 aEdge = TopoDS::Edge(occgeo.emap(aGeomEdgeInd));
2081 param = seg.epgeominfo[j].dist;
2084 else // middle point
2086 param = param2 * 0.5;
2088 if (!aEdge.IsNull() && nodeVec_ACCESS(pind)->getshapeId() < 1)
2090 meshDS->SetNodeOnEdge(nodeVec_ACCESS(pind), aEdge, param);
2095 SMDS_MeshEdge* edge = 0;
2096 if (nbp == 2) // second order ?
2098 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1])))
2100 if ( quadHelper ) // final mesh must be quadratic
2101 edge = quadHelper->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
2103 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
2107 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
2108 nodeVec_ACCESS(pinds[2])))
2110 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
2111 nodeVec_ACCESS(pinds[2]));
2115 if ( comment.empty() ) comment << "Cannot create a mesh edge";
2116 MESSAGE("Cannot create a mesh edge");
2117 nbSeg = nbFac = nbVol = 0;
2120 if ( !aEdge.IsNull() && edge->getshapeId() < 1 )
2121 meshDS->SetMeshElementOnShape(edge, aEdge);
2123 else if ( comment.empty() )
2125 comment << "Invalid netgen segment #" << i;
2129 // ----------------------------------------
2130 // Create mesh faces along geometric faces
2131 // ----------------------------------------
2133 int nbInitFac = initState._nbFaces;
2134 int quadFaceID = ngMesh.GetNFD() + 1;
2135 if ( nbInitFac < nbFac )
2136 // add a faces descriptor to exclude qudrangle elements generated by NETGEN
2137 // from computation of 3D mesh
2138 ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(quadFaceID, /*solid1=*/0, /*solid2=*/0, 0));
2140 vector<const SMDS_MeshNode*> nodes;
2141 for (i = nbInitFac+1; i <= nbFac; ++i )
2143 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
2144 int aGeomFaceInd = elem.GetIndex();
2146 if (aGeomFaceInd > 0 && aGeomFaceInd <= occgeo.fmap.Extent())
2147 aFace = TopoDS::Face(occgeo.fmap(aGeomFaceInd));
2149 for (int j=1; j <= elem.GetNP(); ++j)
2151 int pind = elem.PNum(j);
2152 if ( pind < 1 || pind >= (int) nodeVec.size() )
2154 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind))
2156 nodes.push_back( node );
2157 if (!aFace.IsNull() && node->getshapeId() < 1)
2159 const netgen::PointGeomInfo& pgi = elem.GeomInfoPi(j);
2160 meshDS->SetNodeOnFace(node, aFace, pgi.u, pgi.v);
2164 if ((int) nodes.size() != elem.GetNP() )
2166 if ( comment.empty() )
2167 comment << "Invalid netgen 2d element #" << i;
2168 continue; // bad node ids
2170 SMDS_MeshFace* face = NULL;
2171 switch (elem.GetType())
2174 if ( quadHelper ) // final mesh must be quadratic
2175 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2]);
2177 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2]);
2180 if ( quadHelper ) // final mesh must be quadratic
2181 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2183 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2184 // exclude qudrangle elements from computation of 3D mesh
2185 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2188 nodes[5] = mediumNode( nodes[0],nodes[1],nodes[5], quadHelper );
2189 nodes[3] = mediumNode( nodes[1],nodes[2],nodes[3], quadHelper );
2190 nodes[4] = mediumNode( nodes[2],nodes[0],nodes[4], quadHelper );
2191 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[5],nodes[3],nodes[4]);
2194 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2195 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2196 nodes[5] = mediumNode( nodes[2],nodes[3],nodes[5], quadHelper );
2197 nodes[6] = mediumNode( nodes[3],nodes[0],nodes[6], quadHelper );
2198 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3],
2199 nodes[4],nodes[7],nodes[5],nodes[6]);
2200 // exclude qudrangle elements from computation of 3D mesh
2201 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2204 MESSAGE("NETGEN created a face of unexpected type, ignoring");
2209 if ( comment.empty() ) comment << "Cannot create a mesh face";
2210 MESSAGE("Cannot create a mesh face");
2211 nbSeg = nbFac = nbVol = 0;
2214 if (!aFace.IsNull())
2215 meshDS->SetMeshElementOnShape(face, aFace);
2218 // ------------------
2219 // Create tetrahedra
2220 // ------------------
2222 for (i = 1; i <= nbVol; ++i)
2224 const netgen::Element& elem = ngMesh.VolumeElement(i);
2225 int aSolidInd = elem.GetIndex();
2226 TopoDS_Solid aSolid;
2227 if (aSolidInd > 0 && aSolidInd <= occgeo.somap.Extent())
2228 aSolid = TopoDS::Solid(occgeo.somap(aSolidInd));
2230 for (int j=1; j <= elem.GetNP(); ++j)
2232 int pind = elem.PNum(j);
2233 if ( pind < 1 || pind >= (int)nodeVec.size() )
2235 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind) )
2237 nodes.push_back(node);
2238 if ( !aSolid.IsNull() && node->getshapeId() < 1 )
2239 meshDS->SetNodeInVolume(node, aSolid);
2242 if ((int) nodes.size() != elem.GetNP() )
2244 if ( comment.empty() )
2245 comment << "Invalid netgen 3d element #" << i;
2248 SMDS_MeshVolume* vol = NULL;
2249 switch (elem.GetType())
2252 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3]);
2255 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2256 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2257 nodes[5] = mediumNode( nodes[2],nodes[0],nodes[5], quadHelper );
2258 nodes[6] = mediumNode( nodes[0],nodes[3],nodes[6], quadHelper );
2259 nodes[8] = mediumNode( nodes[1],nodes[3],nodes[8], quadHelper );
2260 nodes[9] = mediumNode( nodes[2],nodes[3],nodes[9], quadHelper );
2261 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3],
2262 nodes[4],nodes[7],nodes[5],nodes[6],nodes[8],nodes[9]);
2265 MESSAGE("NETGEN created a volume of unexpected type, ignoring");
2270 if ( comment.empty() ) comment << "Cannot create a mesh volume";
2271 MESSAGE("Cannot create a mesh volume");
2272 nbSeg = nbFac = nbVol = 0;
2275 if (!aSolid.IsNull())
2276 meshDS->SetMeshElementOnShape(vol, aSolid);
2278 return comment.empty() ? 0 : 1;
2283 //================================================================================
2285 * \brief Restrict size of elements on the given edge
2287 //================================================================================
2289 void setLocalSize(const TopoDS_Edge& edge,
2293 if ( size <= std::numeric_limits<double>::min() )
2295 Standard_Real u1, u2;
2296 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, u1, u2);
2297 if ( curve.IsNull() )
2299 TopoDS_Iterator vIt( edge );
2300 if ( !vIt.More() ) return;
2301 gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex( vIt.Value() ));
2302 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
2306 const int nb = (int)( 1.5 * SMESH_Algo::EdgeLength( edge ) / size );
2307 Standard_Real delta = (u2-u1)/nb;
2308 for(int i=0; i<nb; i++)
2310 Standard_Real u = u1 + delta*i;
2311 gp_Pnt p = curve->Value(u);
2312 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
2313 netgen::Point3d pi(p.X(), p.Y(), p.Z());
2314 double resultSize = mesh.GetH(pi);
2315 if ( resultSize - size > 0.1*size )
2316 // netgen does restriction iff oldH/newH > 1.2 (localh.cpp:136)
2317 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), resultSize/1.201 );
2322 //================================================================================
2324 * \brief Convert error into text
2326 //================================================================================
2328 std::string text(int err)
2333 SMESH_Comment("Error in netgen::OCCGenerateMesh() at ") << netgen::multithread.task;
2336 //================================================================================
2338 * \brief Convert exception into text
2340 //================================================================================
2342 std::string text(Standard_Failure& ex)
2344 SMESH_Comment str("Exception in netgen::OCCGenerateMesh()");
2345 str << " at " << netgen::multithread.task
2346 << ": " << ex.DynamicType()->Name();
2347 if ( ex.GetMessageString() && strlen( ex.GetMessageString() ))
2348 str << ": " << ex.GetMessageString();
2351 //================================================================================
2353 * \brief Convert exception into text
2355 //================================================================================
2357 std::string text(netgen::NgException& ex)
2359 SMESH_Comment str("NgException");
2360 if ( strlen( netgen::multithread.task ) > 0 )
2361 str << " at " << netgen::multithread.task;
2362 str << ": " << ex.What();
2366 //================================================================================
2368 * \brief Looks for triangles lying on a SOLID
2370 //================================================================================
2372 bool hasBadElemOnSolid( const list<const SMDS_MeshElement*>& elems,
2373 SMESH_subMesh* solidSM )
2375 TopTools_IndexedMapOfShape solidSubs;
2376 TopExp::MapShapes( solidSM->GetSubShape(), solidSubs );
2377 SMESHDS_Mesh* mesh = solidSM->GetFather()->GetMeshDS();
2379 list<const SMDS_MeshElement*>::const_iterator e = elems.begin();
2380 for ( ; e != elems.end(); ++e )
2382 const SMDS_MeshElement* elem = *e;
2383 // if ( elem->GetType() != SMDSAbs_Face ) -- 23047
2385 int nbNodesOnSolid = 0, nbNodes = elem->NbNodes();
2386 SMDS_NodeIteratorPtr nIt = elem->nodeIterator();
2387 while ( nIt->more() )
2389 const SMDS_MeshNode* n = nIt->next();
2390 const TopoDS_Shape& s = mesh->IndexToShape( n->getshapeId() );
2391 nbNodesOnSolid += ( !s.IsNull() && solidSubs.Contains( s ));
2392 if ( nbNodesOnSolid > 2 ||
2393 nbNodesOnSolid == nbNodes)
2400 const double edgeMeshingTime = 0.001;
2401 const double faceMeshingTime = 0.019;
2402 const double edgeFaceMeshingTime = edgeMeshingTime + faceMeshingTime;
2403 const double faceOptimizTime = 0.06;
2404 const double voluMeshingTime = 0.15;
2405 const double volOptimizeTime = 0.77;
2408 //=============================================================================
2410 * Here we are going to use the NETGEN mesher
2412 //=============================================================================
2414 bool NETGENPlugin_Mesher::Compute()
2416 NETGENPlugin_NetgenLibWrapper ngLib;
2418 netgen::MeshingParameters& mparams = netgen::mparam;
2419 MESSAGE("Compute with:\n"
2420 " max size = " << mparams.maxh << "\n"
2421 " segments per edge = " << mparams.segmentsperedge);
2423 " growth rate = " << mparams.grading << "\n"
2424 " elements per radius = " << mparams.curvaturesafety << "\n"
2425 " second order = " << mparams.secondorder << "\n"
2426 " quad allowed = " << mparams.quad << "\n"
2427 " surface curvature = " << mparams.uselocalh << "\n"
2428 " fuse edges = " << netgen::merge_solids);
2430 SMESH_ComputeErrorPtr error = SMESH_ComputeError::New();
2431 SMESH_MesherHelper quadHelper( *_mesh );
2432 quadHelper.SetIsQuadratic( mparams.secondorder );
2434 static string debugFile = "/tmp/ngMesh.py"; /* to call toPython( _ngMesh, debugFile )
2435 while debugging netgen */
2436 // -------------------------
2437 // Prepare OCC geometry
2438 // -------------------------
2440 netgen::OCCGeometry occgeo;
2441 list< SMESH_subMesh* > meshedSM[3]; // for 0-2 dimensions
2442 NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
2443 PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
2446 _totalTime = edgeFaceMeshingTime;
2448 _totalTime += faceOptimizTime;
2450 _totalTime += voluMeshingTime + ( _optimize ? volOptimizeTime : 0 );
2451 double doneTime = 0;
2454 _curShapeIndex = -1;
2456 // -------------------------
2457 // Generate the mesh
2458 // -------------------------
2461 NETGENPlugin_ngMeshInfo initState; // it remembers size of ng mesh equal to size of Smesh
2463 SMESH_Comment comment;
2466 // vector of nodes in which node index == netgen ID
2467 vector< const SMDS_MeshNode* > nodeVec;
2475 // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
2476 mparams.uselocalh = false;
2477 mparams.grading = 0.8; // not limitited size growth
2479 if ( _simpleHyp->GetNumberOfSegments() )
2481 mparams.maxh = occgeo.boundingbox.Diam();
2484 mparams.maxh = _simpleHyp->GetLocalLength();
2487 if ( mparams.maxh == 0.0 )
2488 mparams.maxh = occgeo.boundingbox.Diam();
2489 if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
2490 mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
2492 // Local size on faces
2493 occgeo.face_maxh = mparams.maxh;
2495 // Let netgen create _ngMesh and calculate element size on not meshed shapes
2499 int startWith = netgen::MESHCONST_ANALYSE;
2500 int endWith = netgen::MESHCONST_ANALYSE;
2505 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2507 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2509 if(netgen::multithread.terminate)
2512 comment << text(err);
2514 catch (Standard_Failure& ex)
2516 comment << text(ex);
2518 err = 0; //- MESHCONST_ANALYSE isn't so important step
2521 ngLib.setMesh(( Ng_Mesh*) _ngMesh );
2523 _ngMesh->ClearFaceDescriptors(); // we make descriptors our-self
2527 // Pass 1D simple parameters to NETGEN
2528 // --------------------------------
2529 int nbSeg = _simpleHyp->GetNumberOfSegments();
2530 double segSize = _simpleHyp->GetLocalLength();
2531 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
2533 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
2535 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
2536 setLocalSize( e, segSize, *_ngMesh );
2539 else // if ( ! _simpleHyp )
2541 // Local size on vertices and edges
2542 // --------------------------------
2543 for(std::map<int,double>::const_iterator it=EdgeId2LocalSize.begin(); it!=EdgeId2LocalSize.end(); it++)
2545 int key = (*it).first;
2546 double hi = (*it).second;
2547 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2548 const TopoDS_Edge& e = TopoDS::Edge(shape);
2549 setLocalSize( e, hi, *_ngMesh );
2551 for(std::map<int,double>::const_iterator it=VertexId2LocalSize.begin(); it!=VertexId2LocalSize.end(); it++)
2553 int key = (*it).first;
2554 double hi = (*it).second;
2555 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2556 const TopoDS_Vertex& v = TopoDS::Vertex(shape);
2557 gp_Pnt p = BRep_Tool::Pnt(v);
2558 NETGENPlugin_Mesher::RestrictLocalSize( *_ngMesh, p.XYZ(), hi );
2560 for(map<int,double>::const_iterator it=FaceId2LocalSize.begin();
2561 it!=FaceId2LocalSize.end(); it++)
2563 int key = (*it).first;
2564 double val = (*it).second;
2565 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2566 int faceNgID = occgeo.fmap.FindIndex(shape);
2567 occgeo.SetFaceMaxH(faceNgID, val);
2568 for ( TopExp_Explorer edgeExp( shape, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
2569 setLocalSize( TopoDS::Edge( edgeExp.Current() ), val, *_ngMesh );
2573 // Precompute internal edges (issue 0020676) in order to
2574 // add mesh on them correctly (twice) to netgen mesh
2575 if ( !err && internals.hasInternalEdges() )
2577 // load internal shapes into OCCGeometry
2578 netgen::OCCGeometry intOccgeo;
2579 internals.getInternalEdges( intOccgeo.fmap, intOccgeo.emap, intOccgeo.vmap, meshedSM );
2580 intOccgeo.boundingbox = occgeo.boundingbox;
2581 intOccgeo.shape = occgeo.shape;
2582 intOccgeo.face_maxh.SetSize(intOccgeo.fmap.Extent());
2583 intOccgeo.face_maxh = netgen::mparam.maxh;
2584 netgen::Mesh *tmpNgMesh = NULL;
2588 // compute local H on internal shapes in the main mesh
2589 //OCCSetLocalMeshSize(intOccgeo, *_ngMesh); it deletes _ngMesh->localH
2591 // let netgen create a temporary mesh
2593 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2595 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2597 if(netgen::multithread.terminate)
2600 // copy LocalH from the main to temporary mesh
2601 initState.transferLocalH( _ngMesh, tmpNgMesh );
2603 // compute mesh on internal edges
2604 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2606 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2608 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2610 comment << text(err);
2612 catch (Standard_Failure& ex)
2614 comment << text(ex);
2617 initState.restoreLocalH( tmpNgMesh );
2619 // fill SMESH by netgen mesh
2620 vector< const SMDS_MeshNode* > tmpNodeVec;
2621 FillSMesh( intOccgeo, *tmpNgMesh, initState, *_mesh, tmpNodeVec, comment );
2622 err = ( err || !comment.empty() );
2624 nglib::Ng_DeleteMesh((nglib::Ng_Mesh*)tmpNgMesh);
2627 // Fill _ngMesh with nodes and segments of computed submeshes
2630 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_0D ]) &&
2631 FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_1D ], &quadHelper));
2633 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2638 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2643 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2645 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2647 if(netgen::multithread.terminate)
2650 comment << text(err);
2652 catch (Standard_Failure& ex)
2654 comment << text(ex);
2659 _ticTime = ( doneTime += edgeMeshingTime ) / _totalTime / _progressTic;
2661 mparams.uselocalh = true; // restore as it is used at surface optimization
2663 // ---------------------
2664 // compute surface mesh
2665 // ---------------------
2668 // Pass 2D simple parameters to NETGEN
2670 if ( double area = _simpleHyp->GetMaxElementArea() ) {
2672 mparams.maxh = sqrt(2. * area/sqrt(3.0));
2673 mparams.grading = 0.4; // moderate size growth
2676 // length from edges
2677 if ( _ngMesh->GetNSeg() ) {
2678 double edgeLength = 0;
2679 TopTools_MapOfShape visitedEdges;
2680 for ( TopExp_Explorer exp( _shape, TopAbs_EDGE ); exp.More(); exp.Next() )
2681 if( visitedEdges.Add(exp.Current()) )
2682 edgeLength += SMESH_Algo::EdgeLength( TopoDS::Edge( exp.Current() ));
2683 // we have to multiply length by 2 since for each TopoDS_Edge there
2684 // are double set of NETGEN edges, in other words, we have to
2685 // divide _ngMesh->GetNSeg() by 2.
2686 mparams.maxh = 2*edgeLength / _ngMesh->GetNSeg();
2689 mparams.maxh = 1000;
2691 mparams.grading = 0.2; // slow size growth
2693 mparams.quad = _simpleHyp->GetAllowQuadrangles();
2694 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2695 _ngMesh->SetGlobalH (mparams.maxh);
2696 netgen::Box<3> bb = occgeo.GetBoundingBox();
2697 bb.Increase (bb.Diam()/20);
2698 _ngMesh->SetLocalH (bb.PMin(), bb.PMax(), mparams.grading);
2701 // Care of vertices internal in faces (issue 0020676)
2702 if ( internals.hasInternalVertexInFace() )
2704 // store computed segments in SMESH in order not to create SMESH
2705 // edges for ng segments added by AddIntVerticesInFaces()
2706 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2707 // add segments to faces with internal vertices
2708 AddIntVerticesInFaces( occgeo, *_ngMesh, nodeVec, internals );
2709 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2712 // Build viscous layers
2713 if ( _isViscousLayers2D )
2715 if ( !internals.hasInternalVertexInFace() ) {
2716 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2717 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2719 SMESH_ProxyMesh::Ptr viscousMesh;
2720 SMESH_MesherHelper helper( *_mesh );
2721 for ( int faceID = 1; faceID <= occgeo.fmap.Extent(); ++faceID )
2723 const TopoDS_Face& F = TopoDS::Face( occgeo.fmap( faceID ));
2724 viscousMesh = StdMeshers_ViscousLayers2D::Compute( *_mesh, F );
2727 // exclude from computation ng segments built on EDGEs of F
2728 for (int i = 1; i <= _ngMesh->GetNSeg(); i++)
2730 netgen::Segment & seg = _ngMesh->LineSegment(i);
2731 if (seg.si == faceID)
2734 // add new segments to _ngMesh instead of excluded ones
2735 helper.SetSubShape( F );
2737 StdMeshers_FaceSide::GetFaceWires( F, *_mesh, /*skipMediumNodes=*/true,
2738 error, viscousMesh );
2739 error = AddSegmentsToMesh( *_ngMesh, occgeo, wires, helper, nodeVec );
2741 if ( !error ) error = SMESH_ComputeError::New();
2743 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2746 // Let netgen compute 2D mesh
2747 startWith = netgen::MESHCONST_MESHSURFACE;
2748 endWith = _optimize ? netgen::MESHCONST_OPTSURFACE : netgen::MESHCONST_MESHSURFACE;
2753 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2755 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2757 if(netgen::multithread.terminate)
2760 comment << text (err);
2762 catch (Standard_Failure& ex)
2764 comment << text(ex);
2765 //err = 1; -- try to make volumes anyway
2767 catch (netgen::NgException exc)
2769 comment << text(exc);
2770 //err = 1; -- try to make volumes anyway
2775 doneTime += faceMeshingTime + ( _optimize ? faceOptimizTime : 0 );
2776 _ticTime = doneTime / _totalTime / _progressTic;
2778 // ---------------------
2779 // generate volume mesh
2780 // ---------------------
2781 // Fill _ngMesh with nodes and faces of computed 2D submeshes
2782 if ( !err && _isVolume && ( !meshedSM[ MeshDim_2D ].empty() || mparams.quad ))
2784 // load SMESH with computed segments and faces
2785 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2787 // compute pyramids on quadrangles
2788 SMESH_ProxyMesh::Ptr proxyMesh;
2789 if ( _mesh->NbQuadrangles() > 0 )
2790 for ( int iS = 1; iS <= occgeo.somap.Extent(); ++iS )
2792 StdMeshers_QuadToTriaAdaptor* Adaptor = new StdMeshers_QuadToTriaAdaptor;
2793 proxyMesh.reset( Adaptor );
2795 int nbPyrams = _mesh->NbPyramids();
2796 Adaptor->Compute( *_mesh, occgeo.somap(iS) );
2797 if ( nbPyrams != _mesh->NbPyramids() )
2799 list< SMESH_subMesh* > quadFaceSM;
2800 for (TopExp_Explorer face(occgeo.somap(iS), TopAbs_FACE); face.More(); face.Next())
2801 if ( Adaptor->GetProxySubMesh( face.Current() ))
2803 quadFaceSM.push_back( _mesh->GetSubMesh( face.Current() ));
2804 meshedSM[ MeshDim_2D ].remove( quadFaceSM.back() );
2806 FillNgMesh(occgeo, *_ngMesh, nodeVec, quadFaceSM, &quadHelper, proxyMesh);
2809 // fill _ngMesh with faces of sub-meshes
2810 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_2D ], &quadHelper));
2811 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2812 //toPython( _ngMesh, "/tmp/ngPython.py");
2814 if (!err && _isVolume)
2816 // Pass 3D simple parameters to NETGEN
2817 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
2818 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
2820 if ( double vol = simple3d->GetMaxElementVolume() ) {
2822 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
2823 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2826 // length from faces
2827 mparams.maxh = _ngMesh->AverageH();
2829 _ngMesh->SetGlobalH (mparams.maxh);
2830 mparams.grading = 0.4;
2832 _ngMesh->CalcLocalH(mparams.grading);
2834 _ngMesh->CalcLocalH();
2837 // Care of vertices internal in solids and internal faces (issue 0020676)
2838 if ( internals.hasInternalVertexInSolid() || internals.hasInternalFaces() )
2840 // store computed faces in SMESH in order not to create SMESH
2841 // faces for ng faces added here
2842 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2843 // add ng faces to solids with internal vertices
2844 AddIntVerticesInSolids( occgeo, *_ngMesh, nodeVec, internals );
2845 // duplicate mesh faces on internal faces
2846 FixIntFaces( occgeo, *_ngMesh, internals );
2847 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2849 // Let netgen compute 3D mesh
2850 startWith = endWith = netgen::MESHCONST_MESHVOLUME;
2855 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2857 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2859 if(netgen::multithread.terminate)
2862 if ( comment.empty() ) // do not overwrite a previos error
2863 comment << text(err);
2865 catch (Standard_Failure& ex)
2867 if ( comment.empty() ) // do not overwrite a previos error
2868 comment << text(ex);
2871 catch (netgen::NgException exc)
2873 if ( comment.empty() ) // do not overwrite a previos error
2874 comment << text(exc);
2877 _ticTime = ( doneTime += voluMeshingTime ) / _totalTime / _progressTic;
2879 // Let netgen optimize 3D mesh
2880 if ( !err && _optimize )
2882 startWith = endWith = netgen::MESHCONST_OPTVOLUME;
2887 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2889 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2891 if(netgen::multithread.terminate)
2894 if ( comment.empty() ) // do not overwrite a previos error
2895 comment << text(err);
2897 catch (Standard_Failure& ex)
2899 if ( comment.empty() ) // do not overwrite a previos error
2900 comment << text(ex);
2902 catch (netgen::NgException exc)
2904 if ( comment.empty() ) // do not overwrite a previos error
2905 comment << text(exc);
2909 if (!err && mparams.secondorder > 0)
2914 if ( !meshedSM[ MeshDim_1D ].empty() )
2916 // remove segments not attached to geometry (IPAL0052479)
2917 for (int i = 1; i <= _ngMesh->GetNSeg(); ++i)
2919 const netgen::Segment & seg = _ngMesh->LineSegment (i);
2920 if ( seg.epgeominfo[ 0 ].edgenr == 0 )
2921 _ngMesh->DeleteSegment( i );
2923 _ngMesh->Compress();
2925 // convert to quadratic
2926 netgen::OCCRefinementSurfaces ref (occgeo);
2927 ref.MakeSecondOrder (*_ngMesh);
2929 // care of elements already loaded to SMESH
2930 // if ( initState._nbSegments > 0 )
2931 // makeQuadratic( occgeo.emap, _mesh );
2932 // if ( initState._nbFaces > 0 )
2933 // makeQuadratic( occgeo.fmap, _mesh );
2935 catch (Standard_Failure& ex)
2937 if ( comment.empty() ) // do not overwrite a previos error
2938 comment << "Exception in netgen at passing to 2nd order ";
2940 catch (netgen::NgException exc)
2942 if ( comment.empty() ) // do not overwrite a previos error
2943 comment << exc.What();
2948 _ticTime = 0.98 / _progressTic;
2950 int nbNod = _ngMesh->GetNP();
2951 int nbSeg = _ngMesh->GetNSeg();
2952 int nbFac = _ngMesh->GetNSE();
2953 int nbVol = _ngMesh->GetNE();
2954 bool isOK = ( !err && (_isVolume ? (nbVol > 0) : (nbFac > 0)) );
2956 MESSAGE((err ? "Mesh Generation failure" : "End of Mesh Generation") <<
2957 ", nb nodes: " << nbNod <<
2958 ", nb segments: " << nbSeg <<
2959 ", nb faces: " << nbFac <<
2960 ", nb volumes: " << nbVol);
2962 // Feed back the SMESHDS with the generated Nodes and Elements
2963 if ( true /*isOK*/ ) // get whatever built
2965 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2967 if ( quadHelper.GetIsQuadratic() ) // remove free nodes
2968 for ( size_t i = 0; i < nodeVec.size(); ++i )
2969 if ( nodeVec[i] && nodeVec[i]->NbInverseElements() == 0 )
2970 _mesh->GetMeshDS()->RemoveFreeNode( nodeVec[i], 0, /*fromGroups=*/false );
2972 SMESH_ComputeErrorPtr readErr = ReadErrors(nodeVec);
2973 if ( readErr && !readErr->myBadElements.empty() )
2976 if ( !comment.empty() && !readErr->myComment.empty() ) comment += "\n";
2977 comment += readErr->myComment;
2979 if ( error->IsOK() && ( !isOK || comment.size() > 0 ))
2980 error->myName = COMPERR_ALGO_FAILED;
2981 if ( !comment.empty() )
2982 error->myComment = comment;
2984 // SetIsAlwaysComputed( true ) to empty sub-meshes, which
2985 // appear if the geometry contains coincident sub-shape due
2986 // to bool merge_solids = 1; in netgen/libsrc/occ/occgenmesh.cpp
2987 const int nbMaps = 2;
2988 const TopTools_IndexedMapOfShape* geoMaps[nbMaps] =
2989 { & occgeo.vmap, & occgeo.emap/*, & occgeo.fmap*/ };
2990 for ( int iMap = 0; iMap < nbMaps; ++iMap )
2991 for (int i = 1; i <= geoMaps[iMap]->Extent(); i++)
2992 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( geoMaps[iMap]->FindKey(i)))
2993 if ( !sm->IsMeshComputed() )
2994 sm->SetIsAlwaysComputed( true );
2996 // set bad compute error to subshapes of all failed sub-shapes
2997 if ( !error->IsOK() )
2999 bool pb2D = false, pb3D = false;
3000 for (int i = 1; i <= occgeo.fmap.Extent(); i++) {
3001 int status = occgeo.facemeshstatus[i-1];
3002 if (status == netgen::FACE_MESHED_OK ) continue;
3003 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.fmap( i ))) {
3004 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
3005 if ( !smError || smError->IsOK() ) {
3006 if ( status == netgen::FACE_FAILED )
3007 smError.reset( new SMESH_ComputeError( *error ));
3009 smError.reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED, "Ignored" ));
3010 if ( SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
3011 smError->myName = COMPERR_WARNING;
3013 pb2D = pb2D || smError->IsKO();
3016 if ( !pb2D ) // all faces are OK
3017 for (int i = 1; i <= occgeo.somap.Extent(); i++)
3018 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.somap( i )))
3020 bool smComputed = nbVol && !sm->IsEmpty();
3021 if ( smComputed && internals.hasInternalVertexInSolid( sm->GetId() ))
3023 int nbIntV = internals.getSolidsWithVertices().find( sm->GetId() )->second.size();
3024 SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
3025 smComputed = ( smDS->NbElements() > 0 || smDS->NbNodes() > nbIntV );
3027 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
3028 if ( !smComputed && ( !smError || smError->IsOK() ))
3030 smError.reset( new SMESH_ComputeError( *error ));
3031 if ( nbVol && SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
3033 smError->myName = COMPERR_WARNING;
3035 else if ( !smError->myBadElements.empty() ) // bad surface mesh
3037 if ( !hasBadElemOnSolid( smError->myBadElements, sm ))
3041 pb3D = pb3D || ( smError && smError->IsKO() );
3043 if ( !pb2D && !pb3D )
3044 err = 0; // no fatal errors, only warnings
3047 ngLib._isComputeOk = !err;
3052 //=============================================================================
3056 //=============================================================================
3057 bool NETGENPlugin_Mesher::Evaluate(MapShapeNbElems& aResMap)
3059 netgen::MeshingParameters& mparams = netgen::mparam;
3062 // -------------------------
3063 // Prepare OCC geometry
3064 // -------------------------
3065 netgen::OCCGeometry occgeo;
3066 list< SMESH_subMesh* > meshedSM[4]; // for 0-3 dimensions
3067 NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
3068 PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
3070 bool tooManyElems = false;
3071 const int hugeNb = std::numeric_limits<int>::max() / 100;
3076 // pass 1D simple parameters to NETGEN
3079 // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
3080 mparams.uselocalh = false;
3081 mparams.grading = 0.8; // not limitited size growth
3083 if ( _simpleHyp->GetNumberOfSegments() )
3085 mparams.maxh = occgeo.boundingbox.Diam();
3088 mparams.maxh = _simpleHyp->GetLocalLength();
3091 if ( mparams.maxh == 0.0 )
3092 mparams.maxh = occgeo.boundingbox.Diam();
3093 if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
3094 mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
3096 // let netgen create _ngMesh and calculate element size on not meshed shapes
3097 NETGENPlugin_NetgenLibWrapper ngLib;
3098 netgen::Mesh *ngMesh = NULL;
3102 int startWith = netgen::MESHCONST_ANALYSE;
3103 int endWith = netgen::MESHCONST_MESHEDGES;
3105 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, mparams, startWith, endWith);
3107 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, startWith, endWith, optstr);
3110 if(netgen::multithread.terminate)
3113 ngLib.setMesh(( Ng_Mesh*) ngMesh );
3115 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( _shape ))
3116 sm->GetComputeError().reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED ));
3121 // Pass 1D simple parameters to NETGEN
3122 // --------------------------------
3123 int nbSeg = _simpleHyp->GetNumberOfSegments();
3124 double segSize = _simpleHyp->GetLocalLength();
3125 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
3127 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
3129 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
3130 setLocalSize( e, segSize, *ngMesh );
3133 else // if ( ! _simpleHyp )
3135 // Local size on vertices and edges
3136 // --------------------------------
3137 for(std::map<int,double>::const_iterator it=EdgeId2LocalSize.begin(); it!=EdgeId2LocalSize.end(); it++)
3139 int key = (*it).first;
3140 double hi = (*it).second;
3141 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3142 const TopoDS_Edge& e = TopoDS::Edge(shape);
3143 setLocalSize( e, hi, *ngMesh );
3145 for(std::map<int,double>::const_iterator it=VertexId2LocalSize.begin(); it!=VertexId2LocalSize.end(); it++)
3147 int key = (*it).first;
3148 double hi = (*it).second;
3149 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3150 const TopoDS_Vertex& v = TopoDS::Vertex(shape);
3151 gp_Pnt p = BRep_Tool::Pnt(v);
3152 NETGENPlugin_Mesher::RestrictLocalSize( *ngMesh, p.XYZ(), hi );
3154 for(map<int,double>::const_iterator it=FaceId2LocalSize.begin();
3155 it!=FaceId2LocalSize.end(); it++)
3157 int key = (*it).first;
3158 double val = (*it).second;
3159 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3160 int faceNgID = occgeo.fmap.FindIndex(shape);
3161 occgeo.SetFaceMaxH(faceNgID, val);
3162 for ( TopExp_Explorer edgeExp( shape, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
3163 setLocalSize( TopoDS::Edge( edgeExp.Current() ), val, *ngMesh );
3166 // calculate total nb of segments and length of edges
3167 double fullLen = 0.0;
3169 int entity = mparams.secondorder > 0 ? SMDSEntity_Quad_Edge : SMDSEntity_Edge;
3170 TopTools_DataMapOfShapeInteger Edge2NbSeg;
3171 for (TopExp_Explorer exp(_shape, TopAbs_EDGE); exp.More(); exp.Next())
3173 TopoDS_Edge E = TopoDS::Edge( exp.Current() );
3174 if( !Edge2NbSeg.Bind(E,0) )
3177 double aLen = SMESH_Algo::EdgeLength(E);
3180 vector<int>& aVec = aResMap[_mesh->GetSubMesh(E)];
3182 aVec.resize( SMDSEntity_Last, 0);
3184 fullNbSeg += aVec[ entity ];
3187 // store nb of segments computed by Netgen
3188 NCollection_Map<Link> linkMap;
3189 for (int i = 1; i <= ngMesh->GetNSeg(); ++i )
3191 const netgen::Segment& seg = ngMesh->LineSegment(i);
3192 Link link(seg[0], seg[1]);
3193 if ( !linkMap.Add( link )) continue;
3194 int aGeomEdgeInd = seg.epgeominfo[0].edgenr;
3195 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
3197 vector<int>& aVec = aResMap[_mesh->GetSubMesh(occgeo.emap(aGeomEdgeInd))];
3201 // store nb of nodes on edges computed by Netgen
3202 TopTools_DataMapIteratorOfDataMapOfShapeInteger Edge2NbSegIt(Edge2NbSeg);
3203 for (; Edge2NbSegIt.More(); Edge2NbSegIt.Next())
3205 vector<int>& aVec = aResMap[_mesh->GetSubMesh(Edge2NbSegIt.Key())];
3206 if ( aVec[ entity ] > 1 && aVec[ SMDSEntity_Node ] == 0 )
3207 aVec[SMDSEntity_Node] = mparams.secondorder > 0 ? 2*aVec[ entity ]-1 : aVec[ entity ]-1;
3209 fullNbSeg += aVec[ entity ];
3210 Edge2NbSeg( Edge2NbSegIt.Key() ) = aVec[ entity ];
3212 if ( fullNbSeg == 0 )
3219 if ( double area = _simpleHyp->GetMaxElementArea() ) {
3221 mparams.maxh = sqrt(2. * area/sqrt(3.0));
3222 mparams.grading = 0.4; // moderate size growth
3225 // length from edges
3226 mparams.maxh = fullLen/fullNbSeg;
3227 mparams.grading = 0.2; // slow size growth
3230 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3231 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3233 for (TopExp_Explorer exp(_shape, TopAbs_FACE); exp.More(); exp.Next())
3235 TopoDS_Face F = TopoDS::Face( exp.Current() );
3236 SMESH_subMesh *sm = _mesh->GetSubMesh(F);
3238 BRepGProp::SurfaceProperties(F,G);
3239 double anArea = G.Mass();
3240 tooManyElems = tooManyElems || ( anArea/hugeNb > mparams.maxh*mparams.maxh );
3242 if ( !tooManyElems )
3244 TopTools_MapOfShape egdes;
3245 for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next())
3246 if ( egdes.Add( exp1.Current() ))
3247 nb1d += Edge2NbSeg.Find(exp1.Current());
3249 int nbFaces = tooManyElems ? hugeNb : int( 4*anArea / (mparams.maxh*mparams.maxh*sqrt(3.)));
3250 int nbNodes = tooManyElems ? hugeNb : (( nbFaces*3 - (nb1d-1)*2 ) / 6 + 1 );
3252 vector<int> aVec(SMDSEntity_Last, 0);
3253 if( mparams.secondorder > 0 ) {
3254 int nb1d_in = (nbFaces*3 - nb1d) / 2;
3255 aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
3256 aVec[SMDSEntity_Quad_Triangle] = nbFaces;
3259 aVec[SMDSEntity_Node] = Max ( nbNodes, 0 );
3260 aVec[SMDSEntity_Triangle] = nbFaces;
3262 aResMap[sm].swap(aVec);
3269 // pass 3D simple parameters to NETGEN
3270 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
3271 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
3273 if ( double vol = simple3d->GetMaxElementVolume() ) {
3275 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
3276 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3279 // using previous length from faces
3281 mparams.grading = 0.4;
3282 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3285 BRepGProp::VolumeProperties(_shape,G);
3286 double aVolume = G.Mass();
3287 double tetrVol = 0.1179*mparams.maxh*mparams.maxh*mparams.maxh;
3288 tooManyElems = tooManyElems || ( aVolume/hugeNb > tetrVol );
3289 int nbVols = tooManyElems ? hugeNb : int(aVolume/tetrVol);
3290 int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
3291 vector<int> aVec(SMDSEntity_Last, 0 );
3292 if ( tooManyElems ) // avoid FPE
3294 aVec[SMDSEntity_Node] = hugeNb;
3295 aVec[ mparams.secondorder > 0 ? SMDSEntity_Quad_Tetra : SMDSEntity_Tetra] = hugeNb;
3299 if( mparams.secondorder > 0 ) {
3300 aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
3301 aVec[SMDSEntity_Quad_Tetra] = nbVols;
3304 aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
3305 aVec[SMDSEntity_Tetra] = nbVols;
3308 SMESH_subMesh *sm = _mesh->GetSubMesh(_shape);
3309 aResMap[sm].swap(aVec);
3315 double NETGENPlugin_Mesher::GetProgress(const SMESH_Algo* holder,
3316 const int * algoProgressTic,
3317 const double * algoProgress) const
3319 ((int&) _progressTic ) = *algoProgressTic + 1;
3321 if ( !_occgeom ) return 0;
3323 double progress = -1;
3326 if ( _ticTime < 0 && netgen::multithread.task[0] == 'O'/*Optimizing surface*/ )
3328 ((double&) _ticTime ) = edgeFaceMeshingTime / _totalTime / _progressTic;
3330 else if ( !_optimize /*&& _occgeom->fmap.Extent() > 1*/ )
3332 int doneShapeIndex = -1;
3333 while ( doneShapeIndex+1 < _occgeom->facemeshstatus.Size() &&
3334 _occgeom->facemeshstatus[ doneShapeIndex+1 ])
3336 if ( doneShapeIndex+1 != _curShapeIndex )
3338 ((int&) _curShapeIndex) = doneShapeIndex+1;
3339 double doneShapeRate = _curShapeIndex / double( _occgeom->fmap.Extent() );
3340 double doneTime = edgeMeshingTime + doneShapeRate * faceMeshingTime;
3341 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3342 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3343 // << " " << doneTime / _totalTime / _progressTic << endl;
3347 else if ( !_optimize && _occgeom->somap.Extent() > 1 )
3349 int curShapeIndex = _curShapeIndex;
3350 if ( _ngMesh->GetNE() > 0 )
3352 netgen::Element el = (*_ngMesh)[netgen::ElementIndex( _ngMesh->GetNE()-1 )];
3353 curShapeIndex = el.GetIndex();
3355 if ( curShapeIndex != _curShapeIndex )
3357 ((int&) _curShapeIndex) = curShapeIndex;
3358 double doneShapeRate = _curShapeIndex / double( _occgeom->somap.Extent() );
3359 double doneTime = edgeFaceMeshingTime + doneShapeRate * voluMeshingTime;
3360 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3361 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3362 // << " " << doneTime / _totalTime / _progressTic << endl;
3366 progress = Max( *algoProgressTic * _ticTime, *algoProgress );
3369 ((int&) *algoProgressTic )++;
3370 ((double&) *algoProgress) = progress;
3372 //cout << progress << " " << *algoProgressTic << " " << netgen::multithread.task << " "<< _ticTime << endl;
3374 return Min( progress, 0.99 );
3377 //================================================================================
3379 * \brief Remove "test.out" and "problemfaces" files in current directory
3381 //================================================================================
3383 void NETGENPlugin_Mesher::RemoveTmpFiles()
3385 bool rm = SMESH_File("test.out").remove() ;
3387 if (rm && netgen::testout)
3389 delete netgen::testout;
3390 netgen::testout = 0;
3393 SMESH_File("problemfaces").remove();
3394 SMESH_File("occmesh.rep").remove();
3397 //================================================================================
3399 * \brief Read mesh entities preventing successful computation from "test.out" file
3401 //================================================================================
3403 SMESH_ComputeErrorPtr
3404 NETGENPlugin_Mesher::ReadErrors(const vector<const SMDS_MeshNode* >& nodeVec)
3406 SMESH_ComputeErrorPtr err = SMESH_ComputeError::New
3407 (COMPERR_BAD_INPUT_MESH, "Some edges multiple times in surface mesh");
3408 SMESH_File file("test.out");
3410 vector<int> three1(3), three2(3);
3411 const char* badEdgeStr = " multiple times in surface mesh";
3412 const int badEdgeStrLen = strlen( badEdgeStr );
3413 const int nbNodes = nodeVec.size();
3415 while( !file.eof() )
3417 if ( strncmp( file, "Edge ", 5 ) == 0 &&
3418 file.getInts( two ) &&
3419 strncmp( file, badEdgeStr, badEdgeStrLen ) == 0 &&
3420 two[0] < nbNodes && two[1] < nbNodes )
3422 err->myBadElements.push_back( new SMDS_LinearEdge( nodeVec[ two[0]], nodeVec[ two[1]] ));
3423 file += badEdgeStrLen;
3425 else if ( strncmp( file, "Intersecting: ", 14 ) == 0 )
3428 // openelement 18 with open element 126
3432 const char* pos = file;
3433 bool ok = ( strncmp( file, "openelement ", 12 ) == 0 );
3434 ok = ok && file.getInts( two );
3435 ok = ok && file.getInts( three1 );
3436 ok = ok && file.getInts( three2 );
3437 for ( int i = 0; ok && i < 3; ++i )
3438 ok = ( three1[i] < nbNodes && nodeVec[ three1[i]]);
3439 for ( int i = 0; ok && i < 3; ++i )
3440 ok = ( three2[i] < nbNodes && nodeVec[ three2[i]]);
3443 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three1[0]],
3444 nodeVec[ three1[1]],
3445 nodeVec[ three1[2]]));
3446 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three2[0]],
3447 nodeVec[ three2[1]],
3448 nodeVec[ three2[2]]));
3449 err->myComment = "Intersecting triangles";
3463 size_t nbBadElems = err->myBadElements.size();
3470 //================================================================================
3472 * \brief Write a python script creating an equivalent SALOME mesh.
3473 * This is useful to see what mesh is passed as input for the next step of mesh
3474 * generation (of mesh of higher dimension)
3476 //================================================================================
3478 void NETGENPlugin_Mesher::toPython( const netgen::Mesh* ngMesh,
3479 const std::string& pyFile)
3481 ofstream outfile(pyFile.c_str(), ios::out);
3482 if ( !outfile ) return;
3484 outfile << "import SMESH" << endl
3485 << "from salome.smesh import smeshBuilder" << endl
3486 << "smesh = smeshBuilder.New(salome.myStudy)" << endl
3487 << "mesh = smesh.Mesh()" << endl << endl;
3489 using namespace netgen;
3491 for (pi = PointIndex::BASE;
3492 pi < ngMesh->GetNP()+PointIndex::BASE; pi++)
3494 outfile << "mesh.AddNode( ";
3495 outfile << (*ngMesh)[pi](0) << ", ";
3496 outfile << (*ngMesh)[pi](1) << ", ";
3497 outfile << (*ngMesh)[pi](2) << ") ## "<< pi << endl;
3500 int nbDom = ngMesh->GetNDomains();
3501 for ( int i = 0; i < nbDom; ++i )
3502 outfile<< "grp" << i+1 << " = mesh.CreateEmptyGroup( SMESH.FACE, 'domain"<< i+1 << "')"<< endl;
3504 SurfaceElementIndex sei;
3505 for (sei = 0; sei < ngMesh->GetNSE(); sei++)
3507 outfile << "mesh.AddFace([ ";
3508 Element2d sel = (*ngMesh)[sei];
3509 for (int j = 0; j < sel.GetNP(); j++)
3510 outfile << sel[j] << ( j+1 < sel.GetNP() ? ", " : " ])");
3511 if ( sel.IsDeleted() ) outfile << " ## IsDeleted ";
3514 if ((*ngMesh)[sei].GetIndex())
3516 if ( int dom1 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainIn())
3517 outfile << "grp"<< dom1 <<".Add([ " << (int)sei+1 << " ])" << endl;
3518 if ( int dom2 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainOut())
3519 outfile << "grp"<< dom2 <<".Add([ " << (int)sei+1 << " ])" << endl;
3523 for (ElementIndex ei = 0; ei < ngMesh->GetNE(); ei++)
3525 Element el = (*ngMesh)[ei];
3526 outfile << "mesh.AddVolume([ ";
3527 for (int j = 0; j < el.GetNP(); j++)
3528 outfile << el[j] << ( j+1 < el.GetNP() ? ", " : " ])");
3532 for (int i = 1; i <= ngMesh->GetNSeg(); i++)
3534 const Segment & seg = ngMesh->LineSegment (i);
3535 outfile << "mesh.AddEdge([ "
3537 << seg[1] << " ])" << endl;
3539 cout << "Write " << pyFile << endl;
3542 //================================================================================
3544 * \brief Constructor of NETGENPlugin_ngMeshInfo
3546 //================================================================================
3548 NETGENPlugin_ngMeshInfo::NETGENPlugin_ngMeshInfo( netgen::Mesh* ngMesh):
3553 _nbNodes = ngMesh->GetNP();
3554 _nbSegments = ngMesh->GetNSeg();
3555 _nbFaces = ngMesh->GetNSE();
3556 _nbVolumes = ngMesh->GetNE();
3560 _nbNodes = _nbSegments = _nbFaces = _nbVolumes = 0;
3564 //================================================================================
3566 * \brief Copy LocalH member from one netgen mesh to another
3568 //================================================================================
3570 void NETGENPlugin_ngMeshInfo::transferLocalH( netgen::Mesh* fromMesh,
3571 netgen::Mesh* toMesh )
3573 if ( !fromMesh->LocalHFunctionGenerated() ) return;
3574 if ( !toMesh->LocalHFunctionGenerated() )
3576 toMesh->CalcLocalH(netgen::mparam.grading);
3578 toMesh->CalcLocalH();
3581 const size_t size = sizeof( netgen::LocalH );
3582 _copyOfLocalH = new char[ size ];
3583 memcpy( (void*)_copyOfLocalH, (void*)&toMesh->LocalHFunction(), size );
3584 memcpy( (void*)&toMesh->LocalHFunction(), (void*)&fromMesh->LocalHFunction(), size );
3587 //================================================================================
3589 * \brief Restore LocalH member of a netgen mesh
3591 //================================================================================
3593 void NETGENPlugin_ngMeshInfo::restoreLocalH( netgen::Mesh* toMesh )
3595 if ( _copyOfLocalH )
3597 const size_t size = sizeof( netgen::LocalH );
3598 memcpy( (void*)&toMesh->LocalHFunction(), (void*)_copyOfLocalH, size );
3599 delete [] _copyOfLocalH;
3604 //================================================================================
3606 * \brief Find "internal" sub-shapes
3608 //================================================================================
3610 NETGENPlugin_Internals::NETGENPlugin_Internals( SMESH_Mesh& mesh,
3611 const TopoDS_Shape& shape,
3613 : _mesh( mesh ), _is3D( is3D )
3615 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3617 TopExp_Explorer f,e;
3618 for ( f.Init( shape, TopAbs_FACE ); f.More(); f.Next() )
3620 int faceID = meshDS->ShapeToIndex( f.Current() );
3622 // find not computed internal edges
3624 for ( e.Init( f.Current().Oriented(TopAbs_FORWARD), TopAbs_EDGE ); e.More(); e.Next() )
3625 if ( e.Current().Orientation() == TopAbs_INTERNAL )
3627 SMESH_subMesh* eSM = mesh.GetSubMesh( e.Current() );
3628 if ( eSM->IsEmpty() )
3630 _e2face.insert( make_pair( eSM->GetId(), faceID ));
3631 for ( TopoDS_Iterator v(e.Current()); v.More(); v.Next() )
3632 _e2face.insert( make_pair( meshDS->ShapeToIndex( v.Value() ), faceID ));
3636 // find internal vertices in a face
3637 set<int> intVV; // issue 0020850 where same vertex is twice in a face
3638 for ( TopoDS_Iterator fSub( f.Current() ); fSub.More(); fSub.Next())
3639 if ( fSub.Value().ShapeType() == TopAbs_VERTEX )
3641 int vID = meshDS->ShapeToIndex( fSub.Value() );
3642 if ( intVV.insert( vID ).second )
3643 _f2v[ faceID ].push_back( vID );
3648 // find internal faces and their subshapes where nodes are to be doubled
3649 // to make a crack with non-sewed borders
3651 if ( f.Current().Orientation() == TopAbs_INTERNAL )
3653 _intShapes.insert( meshDS->ShapeToIndex( f.Current() ));
3656 list< TopoDS_Shape > edges;
3657 for ( e.Init( f.Current(), TopAbs_EDGE ); e.More(); e.Next())
3658 if ( SMESH_MesherHelper::NbAncestors( e.Current(), mesh, TopAbs_FACE ) > 1 )
3660 _intShapes.insert( meshDS->ShapeToIndex( e.Current() ));
3661 edges.push_back( e.Current() );
3662 // find border faces
3663 PShapeIteratorPtr fIt =
3664 SMESH_MesherHelper::GetAncestors( edges.back(),mesh,TopAbs_FACE );
3665 while ( const TopoDS_Shape* pFace = fIt->next() )
3666 if ( !pFace->IsSame( f.Current() ))
3667 _borderFaces.insert( meshDS->ShapeToIndex( *pFace ));
3670 // we consider vertex internal if it is shared by more than one internal edge
3671 list< TopoDS_Shape >::iterator edge = edges.begin();
3672 for ( ; edge != edges.end(); ++edge )
3673 for ( TopoDS_Iterator v( *edge ); v.More(); v.Next() )
3675 set<int> internalEdges;
3676 PShapeIteratorPtr eIt =
3677 SMESH_MesherHelper::GetAncestors( v.Value(),mesh,TopAbs_EDGE );
3678 while ( const TopoDS_Shape* pEdge = eIt->next() )
3680 int edgeID = meshDS->ShapeToIndex( *pEdge );
3681 if ( isInternalShape( edgeID ))
3682 internalEdges.insert( edgeID );
3684 if ( internalEdges.size() > 1 )
3685 _intShapes.insert( meshDS->ShapeToIndex( v.Value() ));
3689 } // loop on geom faces
3691 // find vertices internal in solids
3694 for ( TopExp_Explorer so(shape, TopAbs_SOLID); so.More(); so.Next())
3696 int soID = meshDS->ShapeToIndex( so.Current() );
3697 for ( TopoDS_Iterator soSub( so.Current() ); soSub.More(); soSub.Next())
3698 if ( soSub.Value().ShapeType() == TopAbs_VERTEX )
3699 _s2v[ soID ].push_back( meshDS->ShapeToIndex( soSub.Value() ));
3704 //================================================================================
3706 * \brief Find mesh faces on non-internal geom faces sharing internal edge
3707 * some nodes of which are to be doubled to make the second border of the "crack"
3709 //================================================================================
3711 void NETGENPlugin_Internals::findBorderElements( TIDSortedElemSet & borderElems )
3713 if ( _intShapes.empty() ) return;
3715 SMESH_Mesh& mesh = const_cast<SMESH_Mesh&>(_mesh);
3716 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3718 // loop on internal geom edges
3719 set<int>::const_iterator intShapeId = _intShapes.begin();
3720 for ( ; intShapeId != _intShapes.end(); ++intShapeId )
3722 const TopoDS_Shape& s = meshDS->IndexToShape( *intShapeId );
3723 if ( s.ShapeType() != TopAbs_EDGE ) continue;
3725 // get internal and non-internal geom faces sharing the internal edge <s>
3727 set<int>::iterator bordFace = _borderFaces.end();
3728 PShapeIteratorPtr faces = SMESH_MesherHelper::GetAncestors( s, _mesh, TopAbs_FACE );
3729 while ( const TopoDS_Shape* pFace = faces->next() )
3731 int faceID = meshDS->ShapeToIndex( *pFace );
3732 if ( isInternalShape( faceID ))
3735 bordFace = _borderFaces.insert( faceID ).first;
3737 if ( bordFace == _borderFaces.end() || !intFace ) continue;
3739 // get all links of mesh faces on internal geom face sharing nodes on edge <s>
3740 set< SMESH_OrientedLink > links; //!< links of faces on internal geom face
3741 list<const SMDS_MeshElement*> suspectFaces[2]; //!< mesh faces on border geom faces
3742 int nbSuspectFaces = 0;
3743 SMESHDS_SubMesh* intFaceSM = meshDS->MeshElements( intFace );
3744 if ( !intFaceSM || intFaceSM->NbElements() == 0 ) continue;
3745 SMESH_subMeshIteratorPtr smIt = mesh.GetSubMesh( s )->getDependsOnIterator(true,true);
3746 while ( smIt->more() )
3748 SMESHDS_SubMesh* sm = smIt->next()->GetSubMeshDS();
3749 if ( !sm ) continue;
3750 SMDS_NodeIteratorPtr nIt = sm->GetNodes();
3751 while ( nIt->more() )
3753 const SMDS_MeshNode* nOnEdge = nIt->next();
3754 SMDS_ElemIteratorPtr fIt = nOnEdge->GetInverseElementIterator(SMDSAbs_Face);
3755 while ( fIt->more() )
3757 const SMDS_MeshElement* f = fIt->next();
3758 const int nbNodes = f->NbCornerNodes();
3759 if ( intFaceSM->Contains( f ))
3761 for ( int i = 0; i < nbNodes; ++i )
3762 links.insert( SMESH_OrientedLink( f->GetNode(i), f->GetNode((i+1)%nbNodes)));
3767 for ( int i = 0; i < nbNodes; ++i )
3768 nbDblNodes += isInternalShape( f->GetNode(i)->getshapeId() );
3770 suspectFaces[ nbDblNodes < 2 ].push_back( f );
3776 // suspectFaces[0] having link with same orientation as mesh faces on
3777 // the internal geom face are <borderElems>. suspectFaces[1] have
3778 // only one node on edge <s>, we decide on them later (at the 2nd loop)
3779 // by links of <borderElems> found at the 1st and 2nd loops
3780 set< SMESH_OrientedLink > borderLinks;
3781 for ( int isPostponed = 0; isPostponed < 2; ++isPostponed )
3783 list<const SMDS_MeshElement*>::iterator fIt = suspectFaces[isPostponed].begin();
3784 for ( int nbF = 0; fIt != suspectFaces[isPostponed].end(); ++fIt, ++nbF )
3786 const SMDS_MeshElement* f = *fIt;
3787 bool isBorder = false, linkFound = false, borderLinkFound = false;
3788 list< SMESH_OrientedLink > faceLinks;
3789 int nbNodes = f->NbCornerNodes();
3790 for ( int i = 0; i < nbNodes; ++i )
3792 SMESH_OrientedLink link( f->GetNode(i), f->GetNode((i+1)%nbNodes));
3793 faceLinks.push_back( link );
3796 set< SMESH_OrientedLink >::iterator foundLink = links.find( link );
3797 if ( foundLink != links.end() )
3800 isBorder = ( foundLink->_reversed == link._reversed );
3801 if ( !isBorder && !isPostponed ) break;
3802 faceLinks.pop_back();
3804 else if ( isPostponed && !borderLinkFound )
3806 foundLink = borderLinks.find( link );
3807 if ( foundLink != borderLinks.end() )
3809 borderLinkFound = true;
3810 isBorder = ( foundLink->_reversed != link._reversed );
3817 borderElems.insert( f );
3818 borderLinks.insert( faceLinks.begin(), faceLinks.end() );
3820 else if ( !linkFound && !borderLinkFound )
3822 suspectFaces[1].push_back( f );
3823 if ( nbF > 2 * nbSuspectFaces )
3824 break; // dead loop protection
3831 //================================================================================
3833 * \brief put internal shapes in maps and fill in submeshes to precompute
3835 //================================================================================
3837 void NETGENPlugin_Internals::getInternalEdges( TopTools_IndexedMapOfShape& fmap,
3838 TopTools_IndexedMapOfShape& emap,
3839 TopTools_IndexedMapOfShape& vmap,
3840 list< SMESH_subMesh* > smToPrecompute[])
3842 if ( !hasInternalEdges() ) return;
3843 map<int,int>::const_iterator ev_face = _e2face.begin();
3844 for ( ; ev_face != _e2face.end(); ++ev_face )
3846 const TopoDS_Shape& ev = _mesh.GetMeshDS()->IndexToShape( ev_face->first );
3847 const TopoDS_Shape& face = _mesh.GetMeshDS()->IndexToShape( ev_face->second );
3849 ( ev.ShapeType() == TopAbs_EDGE ? emap : vmap ).Add( ev );
3851 //cout<<"INTERNAL EDGE or VERTEX "<<ev_face->first<<" on face "<<ev_face->second<<endl;
3853 smToPrecompute[ MeshDim_1D ].push_back( _mesh.GetSubMeshContaining( ev_face->first ));
3857 //================================================================================
3859 * \brief return shapes and submeshes to be meshed and already meshed boundary submeshes
3861 //================================================================================
3863 void NETGENPlugin_Internals::getInternalFaces( TopTools_IndexedMapOfShape& fmap,
3864 TopTools_IndexedMapOfShape& emap,
3865 list< SMESH_subMesh* >& intFaceSM,
3866 list< SMESH_subMesh* >& boundarySM)
3868 if ( !hasInternalFaces() ) return;
3870 // <fmap> and <emap> are for not yet meshed shapes
3871 // <intFaceSM> is for submeshes of faces
3872 // <boundarySM> is for meshed edges and vertices
3877 set<int> shapeIDs ( _intShapes );
3878 if ( !_borderFaces.empty() )
3879 shapeIDs.insert( _borderFaces.begin(), _borderFaces.end() );
3881 set<int>::const_iterator intS = shapeIDs.begin();
3882 for ( ; intS != shapeIDs.end(); ++intS )
3884 SMESH_subMesh* sm = _mesh.GetSubMeshContaining( *intS );
3886 if ( sm->GetSubShape().ShapeType() != TopAbs_FACE ) continue;
3888 intFaceSM.push_back( sm );
3890 // add submeshes of not computed internal faces
3891 if ( !sm->IsEmpty() ) continue;
3893 SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(true,true);
3894 while ( smIt->more() )
3897 const TopoDS_Shape& s = sm->GetSubShape();
3899 if ( sm->IsEmpty() )
3902 switch ( s.ShapeType() ) {
3903 case TopAbs_FACE: fmap.Add ( s ); break;
3904 case TopAbs_EDGE: emap.Add ( s ); break;
3910 if ( s.ShapeType() != TopAbs_FACE )
3911 boundarySM.push_back( sm );
3917 //================================================================================
3919 * \brief Return true if given shape is to be precomputed in order to be correctly
3920 * added to netgen mesh
3922 //================================================================================
3924 bool NETGENPlugin_Internals::isShapeToPrecompute(const TopoDS_Shape& s)
3926 int shapeID = _mesh.GetMeshDS()->ShapeToIndex( s );
3927 switch ( s.ShapeType() ) {
3928 case TopAbs_FACE : break; //return isInternalShape( shapeID ) || isBorderFace( shapeID );
3929 case TopAbs_EDGE : return isInternalEdge( shapeID );
3930 case TopAbs_VERTEX: break;
3936 //================================================================================
3938 * \brief Return SMESH
3940 //================================================================================
3942 SMESH_Mesh& NETGENPlugin_Internals::getMesh() const
3944 return const_cast<SMESH_Mesh&>( _mesh );
3947 //================================================================================
3949 * \brief Initialize netgen library
3951 //================================================================================
3953 NETGENPlugin_NetgenLibWrapper::NETGENPlugin_NetgenLibWrapper()
3957 _isComputeOk = false;
3959 if ( !getenv( "KEEP_NETGEN_OUTPUT" ))
3961 // redirect all netgen output (mycout,myerr,cout) to _outputFileName
3962 _outputFileName = getOutputFileName();
3963 netgen::mycout = new ofstream ( _outputFileName.c_str() );
3964 netgen::myerr = netgen::mycout;
3965 _coutBuffer = std::cout.rdbuf();
3967 cout << "NOTE: netgen output is redirected to file " << _outputFileName << endl;
3969 std::cout.rdbuf( netgen::mycout->rdbuf() );
3973 _ngMesh = Ng_NewMesh();
3976 //================================================================================
3978 * \brief Finish using netgen library
3980 //================================================================================
3982 NETGENPlugin_NetgenLibWrapper::~NETGENPlugin_NetgenLibWrapper()
3984 Ng_DeleteMesh( _ngMesh );
3986 NETGENPlugin_Mesher::RemoveTmpFiles();
3988 std::cout.rdbuf( _coutBuffer );
3995 //================================================================================
3997 * \brief Set netgen mesh to delete at destruction
3999 //================================================================================
4001 void NETGENPlugin_NetgenLibWrapper::setMesh( Ng_Mesh* mesh )
4004 Ng_DeleteMesh( _ngMesh );
4008 //================================================================================
4010 * \brief Return a unique file name
4012 //================================================================================
4014 std::string NETGENPlugin_NetgenLibWrapper::getOutputFileName()
4016 std::string aTmpDir = SALOMEDS_Tool::GetTmpDir();
4018 TCollection_AsciiString aGenericName = (char*)aTmpDir.c_str();
4019 aGenericName += "NETGEN_";
4021 aGenericName += getpid();
4023 aGenericName += _getpid();
4025 aGenericName += "_";
4026 aGenericName += Abs((Standard_Integer)(long) aGenericName.ToCString());
4027 aGenericName += ".out";
4029 return aGenericName.ToCString();
4032 //================================================================================
4034 * \brief Remove file with netgen output
4036 //================================================================================
4038 void NETGENPlugin_NetgenLibWrapper::removeOutputFile()
4040 if ( !_outputFileName.empty() )
4042 if ( netgen::mycout )
4044 delete netgen::mycout;
4048 string tmpDir = SALOMEDS_Tool::GetDirFromPath ( _outputFileName );
4049 string aFileName = SALOMEDS_Tool::GetNameFromPath( _outputFileName ) + ".out";
4050 SALOMEDS::ListOfFileNames_var aFiles = new SALOMEDS::ListOfFileNames;
4052 aFiles[0] = aFileName.c_str();
4054 SALOMEDS_Tool::RemoveTemporaryFiles( tmpDir.c_str(), aFiles.in(), true );