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 mparams.quad = hyp->GetQuadAllowed() ? 1 : 0;
270 _optimize = hyp->GetOptimize();
271 _fineness = hyp->GetFineness();
272 mparams.uselocalh = hyp->GetSurfaceCurvature();
273 netgen::merge_solids = hyp->GetFuseEdges();
276 SMESH_Gen_i* smeshGen_i = SMESH_Gen_i::GetSMESHGen();
277 CORBA::Object_var anObject = smeshGen_i->GetNS()->Resolve("/myStudyManager");
278 SALOMEDS::StudyManager_var aStudyMgr = SALOMEDS::StudyManager::_narrow(anObject);
279 SALOMEDS::Study_var myStudy = aStudyMgr->GetStudyByID(hyp->GetStudyId());
281 const NETGENPlugin_Hypothesis::TLocalSize localSizes = hyp->GetLocalSizesAndEntries();
282 NETGENPlugin_Hypothesis::TLocalSize::const_iterator it = localSizes.begin();
283 for ( ; it != localSizes.end() ; it++)
285 std::string entry = (*it).first;
286 double val = (*it).second;
288 GEOM::GEOM_Object_var aGeomObj;
289 SALOMEDS::SObject_var aSObj = myStudy->FindObjectID( entry.c_str() );
290 if ( !aSObj->_is_nil() ) {
291 CORBA::Object_var obj = aSObj->GetObject();
292 aGeomObj = GEOM::GEOM_Object::_narrow(obj);
295 TopoDS_Shape S = smeshGen_i->GeomObjectToShape( aGeomObj.in() );
296 SetLocalSize(S, val);
301 //=============================================================================
303 * Pass simple parameters to NETGEN
305 //=============================================================================
307 void NETGENPlugin_Mesher::SetParameters(const NETGENPlugin_SimpleHypothesis_2D* hyp)
311 SetDefaultParameters();
314 //=============================================================================
316 * Link - a pair of integer numbers
318 //=============================================================================
322 Link(int _n1, int _n2) : n1(_n1), n2(_n2) {}
323 Link() : n1(0), n2(0) {}
324 bool Contains( int n ) const { return n == n1 || n == n2; }
325 bool IsConnected( const Link& other ) const
327 return (( Contains( other.n1 ) || Contains( other.n2 )) && ( this != &other ));
331 int HashCode(const Link& aLink, int aLimit)
333 return HashCode(aLink.n1 + aLink.n2, aLimit);
336 Standard_Boolean IsEqual(const Link& aLink1, const Link& aLink2)
338 return (( aLink1.n1 == aLink2.n1 && aLink1.n2 == aLink2.n2 ) ||
339 ( aLink1.n1 == aLink2.n2 && aLink1.n2 == aLink2.n1 ));
344 //================================================================================
346 * \brief return id of netgen point corresponding to SMDS node
348 //================================================================================
349 typedef map< const SMDS_MeshNode*, int > TNode2IdMap;
351 int ngNodeId( const SMDS_MeshNode* node,
352 netgen::Mesh& ngMesh,
353 TNode2IdMap& nodeNgIdMap)
355 int newNgId = ngMesh.GetNP() + 1;
357 TNode2IdMap::iterator node_id = nodeNgIdMap.insert( make_pair( node, newNgId )).first;
359 if ( node_id->second == newNgId)
361 #if defined(DUMP_SEGMENTS) || defined(DUMP_TRIANGLES)
362 cout << "Ng " << newNgId << " - " << node;
364 netgen::MeshPoint p( netgen::Point<3> (node->X(), node->Y(), node->Z()) );
365 ngMesh.AddPoint( p );
367 return node_id->second;
370 //================================================================================
372 * \brief Return computed EDGEs connected to the given one
374 //================================================================================
376 list< TopoDS_Edge > getConnectedEdges( const TopoDS_Edge& edge,
377 const TopoDS_Face& face,
378 const set< SMESH_subMesh* > & computedSM,
379 const SMESH_MesherHelper& helper,
380 map< SMESH_subMesh*, set< int > >& addedEdgeSM2Faces)
383 list< TopoDS_Edge > edges;
384 list< int > nbEdgesInWire;
385 /*int nbWires =*/ SMESH_Block::GetOrderedEdges( face, edges, nbEdgesInWire);
387 // find <edge> within <edges>
388 list< TopoDS_Edge >::iterator eItFwd = edges.begin();
389 for ( ; eItFwd != edges.end(); ++eItFwd )
390 if ( edge.IsSame( *eItFwd ))
392 if ( eItFwd == edges.end()) return list< TopoDS_Edge>();
394 if ( eItFwd->Orientation() >= TopAbs_INTERNAL )
396 // connected INTERNAL edges returned from GetOrderedEdges() are wrongly oriented
397 // so treat each INTERNAL edge separately
398 TopoDS_Edge e = *eItFwd;
400 edges.push_back( e );
404 // get all computed EDGEs connected to <edge>
406 list< TopoDS_Edge >::iterator eItBack = eItFwd, ePrev;
407 TopoDS_Vertex vCommon;
408 TopTools_MapOfShape eAdded; // map used not to add a seam edge twice to <edges>
411 // put edges before <edge> to <edges> back
412 while ( edges.begin() != eItFwd )
413 edges.splice( edges.end(), edges, edges.begin() );
417 while ( ++eItFwd != edges.end() )
419 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( *eItFwd );
421 bool connected = TopExp::CommonVertex( *ePrev, *eItFwd, vCommon );
422 bool computed = sm->IsMeshComputed();
423 bool added = addedEdgeSM2Faces[ sm ].count( helper.GetSubShapeID() );
424 bool doubled = !eAdded.Add( *eItFwd );
425 bool orientOK = (( ePrev ->Orientation() < TopAbs_INTERNAL ) ==
426 ( eItFwd->Orientation() < TopAbs_INTERNAL ) );
427 if ( !connected || !computed || !orientOK || added || doubled )
429 // stop advancement; move edges from tail to head
430 while ( edges.back() != *ePrev )
431 edges.splice( edges.begin(), edges, --edges.end() );
437 while ( eItBack != edges.begin() )
441 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( *eItBack );
443 bool connected = TopExp::CommonVertex( *ePrev, *eItBack, vCommon );
444 bool computed = sm->IsMeshComputed();
445 bool added = addedEdgeSM2Faces[ sm ].count( helper.GetSubShapeID() );
446 bool doubled = !eAdded.Add( *eItBack );
447 bool orientOK = (( ePrev ->Orientation() < TopAbs_INTERNAL ) ==
448 ( eItBack->Orientation() < TopAbs_INTERNAL ) );
449 if ( !connected || !computed || !orientOK || added || doubled)
452 edges.erase( edges.begin(), ePrev );
456 if ( edges.front() != edges.back() )
458 // assure that the 1st vertex is meshed
459 TopoDS_Edge eLast = edges.back();
460 while ( !SMESH_Algo::VertexNode( SMESH_MesherHelper::IthVertex( 0, edges.front()), helper.GetMeshDS())
462 edges.front() != eLast )
463 edges.splice( edges.end(), edges, edges.begin() );
468 //================================================================================
470 * \brief Make triangulation of a shape precise enough
472 //================================================================================
474 void updateTriangulation( const TopoDS_Shape& shape )
476 // static set< Poly_Triangulation* > updated;
478 // TopLoc_Location loc;
479 // TopExp_Explorer fExp( shape, TopAbs_FACE );
480 // for ( ; fExp.More(); fExp.Next() )
482 // Handle(Poly_Triangulation) triangulation =
483 // BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
484 // if ( triangulation.IsNull() ||
485 // updated.insert( triangulation.operator->() ).second )
487 // BRepTools::Clean (shape);
490 BRepMesh_IncrementalMesh e(shape, 0.01, true);
492 catch (Standard_Failure)
495 // updated.erase( triangulation.operator->() );
496 // triangulation = BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
497 // updated.insert( triangulation.operator->() );
501 //================================================================================
503 * \brief Returns a medium node either existing in SMESH of created by NETGEN
504 * \param [in] corner1 - corner node 1
505 * \param [in] corner2 - corner node 2
506 * \param [in] defaultMedium - the node created by NETGEN
507 * \param [in] helper - holder of medium nodes existing in SMESH
508 * \return const SMDS_MeshNode* - the result node
510 //================================================================================
512 const SMDS_MeshNode* mediumNode( const SMDS_MeshNode* corner1,
513 const SMDS_MeshNode* corner2,
514 const SMDS_MeshNode* defaultMedium,
515 const SMESH_MesherHelper* helper)
519 TLinkNodeMap::const_iterator l2n =
520 helper->GetTLinkNodeMap().find( SMESH_TLink( corner1, corner2 ));
521 if ( l2n != helper->GetTLinkNodeMap().end() )
522 defaultMedium = l2n->second;
524 return defaultMedium;
527 //================================================================================
529 * \brief Assure that mesh on given shapes is quadratic
531 //================================================================================
533 void makeQuadratic( const TopTools_IndexedMapOfShape& shapes,
536 for ( int i = 1; i <= shapes.Extent(); ++i )
538 SMESHDS_SubMesh* smDS = mesh->GetMeshDS()->MeshElements( shapes(i) );
539 if ( !smDS ) continue;
540 SMDS_ElemIteratorPtr elemIt = smDS->GetElements();
541 if ( !elemIt->more() ) continue;
542 const SMDS_MeshElement* e = elemIt->next();
543 if ( !e || e->IsQuadratic() )
546 TIDSortedElemSet elems;
548 while ( elemIt->more() )
549 elems.insert( elems.end(), elemIt->next() );
551 SMESH_MeshEditor( mesh ).ConvertToQuadratic( /*3d=*/false, elems, /*biQuad=*/false );
557 //================================================================================
559 * \brief Initialize netgen::OCCGeometry with OCCT shape
561 //================================================================================
563 void NETGENPlugin_Mesher::PrepareOCCgeometry(netgen::OCCGeometry& occgeo,
564 const TopoDS_Shape& shape,
566 list< SMESH_subMesh* > * meshedSM,
567 NETGENPlugin_Internals* intern)
569 updateTriangulation( shape );
572 BRepBndLib::Add (shape, bb);
573 double x1,y1,z1,x2,y2,z2;
574 bb.Get (x1,y1,z1,x2,y2,z2);
575 MESSAGE("shape bounding box:\n" <<
576 "(" << x1 << " " << y1 << " " << z1 << ") " <<
577 "(" << x2 << " " << y2 << " " << z2 << ")");
578 netgen::Point<3> p1 = netgen::Point<3> (x1,y1,z1);
579 netgen::Point<3> p2 = netgen::Point<3> (x2,y2,z2);
580 occgeo.boundingbox = netgen::Box<3> (p1,p2);
582 occgeo.shape = shape;
585 // fill maps of shapes of occgeo with not yet meshed subshapes
587 // get root submeshes
588 list< SMESH_subMesh* > rootSM;
589 const int shapeID = mesh.GetMeshDS()->ShapeToIndex( shape );
590 if ( shapeID > 0 ) { // SMESH_subMesh with ID 0 may exist, don't use it!
591 rootSM.push_back( mesh.GetSubMesh( shape ));
594 for ( TopoDS_Iterator it( shape ); it.More(); it.Next() )
595 rootSM.push_back( mesh.GetSubMesh( it.Value() ));
600 // add subshapes of empty submeshes
601 list< SMESH_subMesh* >::iterator rootIt = rootSM.begin(), rootEnd = rootSM.end();
602 for ( ; rootIt != rootEnd; ++rootIt ) {
603 SMESH_subMesh * root = *rootIt;
604 SMESH_subMeshIteratorPtr smIt = root->getDependsOnIterator(/*includeSelf=*/true,
605 /*complexShapeFirst=*/true);
606 // to find a right orientation of subshapes (PAL20462)
607 TopTools_IndexedMapOfShape subShapes;
608 TopExp::MapShapes(root->GetSubShape(), subShapes);
609 while ( smIt->more() )
611 SMESH_subMesh* sm = smIt->next();
612 TopoDS_Shape shape = sm->GetSubShape();
613 totNbFaces += ( shape.ShapeType() == TopAbs_FACE );
614 if ( intern && intern->isShapeToPrecompute( shape ))
616 if ( !meshedSM || sm->IsEmpty() )
618 if ( shape.ShapeType() != TopAbs_VERTEX )
619 shape = subShapes( subShapes.FindIndex( shape ));// shape -> index -> oriented shape
620 if ( shape.Orientation() >= TopAbs_INTERNAL )
621 shape.Orientation( TopAbs_FORWARD ); // isuue 0020676
622 switch ( shape.ShapeType() ) {
623 case TopAbs_FACE : occgeo.fmap.Add( shape ); break;
624 case TopAbs_EDGE : occgeo.emap.Add( shape ); break;
625 case TopAbs_VERTEX: occgeo.vmap.Add( shape ); break;
626 case TopAbs_SOLID :occgeo.somap.Add( shape ); break;
630 // collect submeshes of meshed shapes
633 const int dim = SMESH_Gen::GetShapeDim( shape );
634 meshedSM[ dim ].push_back( sm );
638 occgeo.facemeshstatus.SetSize (totNbFaces);
639 occgeo.facemeshstatus = 0;
640 occgeo.face_maxh_modified.SetSize(totNbFaces);
641 occgeo.face_maxh_modified = 0;
642 occgeo.face_maxh.SetSize(totNbFaces);
643 occgeo.face_maxh = netgen::mparam.maxh;
646 //================================================================================
648 * \brief Return a default min size value suitable for the given geometry.
650 //================================================================================
652 double NETGENPlugin_Mesher::GetDefaultMinSize(const TopoDS_Shape& geom,
653 const double maxSize)
655 updateTriangulation( geom );
659 const int* pi[4] = { &i1, &i2, &i3, &i1 };
662 TopExp_Explorer fExp( geom, TopAbs_FACE );
663 for ( ; fExp.More(); fExp.Next() )
665 Handle(Poly_Triangulation) triangulation =
666 BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
667 if ( triangulation.IsNull() ) continue;
668 const double fTol = BRep_Tool::Tolerance( TopoDS::Face( fExp.Current() ));
669 const TColgp_Array1OfPnt& points = triangulation->Nodes();
670 const Poly_Array1OfTriangle& trias = triangulation->Triangles();
671 for ( int iT = trias.Lower(); iT <= trias.Upper(); ++iT )
673 trias(iT).Get( i1, i2, i3 );
674 for ( int j = 0; j < 3; ++j )
676 double dist2 = points(*pi[j]).SquareDistance( points( *pi[j+1] ));
677 if ( dist2 < minh && fTol*fTol < dist2 )
679 bb.Add( points(*pi[j]));
683 if ( minh > 0.25 * bb.SquareExtent() ) // simple geometry, rough triangulation
685 minh = 1e-3 * sqrt( bb.SquareExtent());
686 //cout << "BND BOX minh = " <<minh << endl;
690 minh = 3 * sqrt( minh ); // triangulation for visualization is rather fine
691 //cout << "TRIANGULATION minh = " <<minh << endl;
693 if ( minh > 0.5 * maxSize )
699 //================================================================================
701 * \brief Restrict size of elements at a given point
703 //================================================================================
705 void NETGENPlugin_Mesher::RestrictLocalSize(netgen::Mesh& ngMesh,
708 const bool overrideMinH)
710 if ( size <= std::numeric_limits<double>::min() )
712 if ( netgen::mparam.minh > size )
716 ngMesh.SetMinimalH( size );
717 netgen::mparam.minh = size;
721 size = netgen::mparam.minh;
724 netgen::Point3d pi(p.X(), p.Y(), p.Z());
725 ngMesh.RestrictLocalH( pi, size );
728 //================================================================================
730 * \brief fill ngMesh with nodes and elements of computed submeshes
732 //================================================================================
734 bool NETGENPlugin_Mesher::FillNgMesh(netgen::OCCGeometry& occgeom,
735 netgen::Mesh& ngMesh,
736 vector<const SMDS_MeshNode*>& nodeVec,
737 const list< SMESH_subMesh* > & meshedSM,
738 SMESH_MesherHelper* quadHelper,
739 SMESH_ProxyMesh::Ptr proxyMesh)
741 TNode2IdMap nodeNgIdMap;
742 for ( size_t i = 1; i < nodeVec.size(); ++i )
743 nodeNgIdMap.insert( make_pair( nodeVec[i], i ));
745 TopTools_MapOfShape visitedShapes;
746 map< SMESH_subMesh*, set< int > > visitedEdgeSM2Faces;
747 set< SMESH_subMesh* > computedSM( meshedSM.begin(), meshedSM.end() );
749 SMESH_MesherHelper helper (*_mesh);
751 int faceNgID = ngMesh.GetNFD();
753 list< SMESH_subMesh* >::const_iterator smIt, smEnd = meshedSM.end();
754 for ( smIt = meshedSM.begin(); smIt != smEnd; ++smIt )
756 SMESH_subMesh* sm = *smIt;
757 if ( !visitedShapes.Add( sm->GetSubShape() ))
760 const SMESHDS_SubMesh * smDS = sm->GetSubMeshDS();
761 if ( !smDS ) continue;
763 switch ( sm->GetSubShape().ShapeType() )
765 case TopAbs_EDGE: { // EDGE
766 // ----------------------
767 TopoDS_Edge geomEdge = TopoDS::Edge( sm->GetSubShape() );
768 if ( geomEdge.Orientation() >= TopAbs_INTERNAL )
769 geomEdge.Orientation( TopAbs_FORWARD ); // issue 0020676
771 // Add ng segments for each not meshed FACE the EDGE bounds
772 PShapeIteratorPtr fIt = helper.GetAncestors( geomEdge, *sm->GetFather(), TopAbs_FACE );
773 while ( const TopoDS_Shape * anc = fIt->next() )
775 faceNgID = occgeom.fmap.FindIndex( *anc );
777 continue; // meshed face
779 int faceSMDSId = helper.GetMeshDS()->ShapeToIndex( *anc );
780 if ( visitedEdgeSM2Faces[ sm ].count( faceSMDSId ))
781 continue; // already treated EDGE
783 TopoDS_Face face = TopoDS::Face( occgeom.fmap( faceNgID ));
784 if ( face.Orientation() >= TopAbs_INTERNAL )
785 face.Orientation( TopAbs_FORWARD ); // issue 0020676
787 // get all meshed EDGEs of the FACE connected to geomEdge (issue 0021140)
788 helper.SetSubShape( face );
789 list< TopoDS_Edge > edges = getConnectedEdges( geomEdge, face, computedSM, helper,
790 visitedEdgeSM2Faces );
792 continue; // wrong ancestor?
794 // find out orientation of <edges> within <face>
795 TopoDS_Edge eNotSeam = edges.front();
796 if ( helper.HasSeam() )
798 list< TopoDS_Edge >::iterator eIt = edges.begin();
799 while ( helper.IsRealSeam( *eIt )) ++eIt;
800 if ( eIt != edges.end() )
803 TopAbs_Orientation fOri = helper.GetSubShapeOri( face, eNotSeam );
804 bool isForwad = ( fOri == eNotSeam.Orientation() || fOri >= TopAbs_INTERNAL );
806 // get all nodes from connected <edges>
807 const bool isQuad = smDS->IsQuadratic();
808 StdMeshers_FaceSide fSide( face, edges, _mesh, isForwad, isQuad );
809 const vector<UVPtStruct>& points = fSide.GetUVPtStruct();
810 if ( points.empty() )
811 return false; // invalid node params?
812 int i, nbSeg = fSide.NbSegments();
814 // remember EDGEs of fSide to treat only once
815 for ( int iE = 0; iE < fSide.NbEdges(); ++iE )
816 visitedEdgeSM2Faces[ helper.GetMesh()->GetSubMesh( fSide.Edge(iE )) ].insert(faceSMDSId);
818 double otherSeamParam = 0;
823 int prevNgId = ngNodeId( points[0].node, ngMesh, nodeNgIdMap );
825 for ( i = 0; i < nbSeg; ++i )
827 const UVPtStruct& p1 = points[ i ];
828 const UVPtStruct& p2 = points[ i+1 ];
830 if ( p1.node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX ) //an EDGE begins
833 if ( helper.IsRealSeam( p1.node->getshapeId() ))
835 TopoDS_Edge e = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
836 isSeam = helper.IsRealSeam( e );
839 otherSeamParam = helper.GetOtherParam( helper.GetPeriodicIndex() & 1 ? p2.u : p2.v );
846 seg[1] = prevNgId = ngNodeId( p2.node, ngMesh, nodeNgIdMap );
847 // node param on curve
848 seg.epgeominfo[ 0 ].dist = p1.param;
849 seg.epgeominfo[ 1 ].dist = p2.param;
851 seg.epgeominfo[ 0 ].u = p1.u;
852 seg.epgeominfo[ 0 ].v = p1.v;
853 seg.epgeominfo[ 1 ].u = p2.u;
854 seg.epgeominfo[ 1 ].v = p2.v;
856 //geomEdge = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
857 //seg.epgeominfo[ 0 ].edgenr = seg.epgeominfo[ 1 ].edgenr = occgeom.emap.FindIndex( geomEdge );
859 //seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
860 seg.si = faceNgID; // = geom.fmap.FindIndex (face);
861 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
862 ngMesh.AddSegment (seg);
864 SMESH_TNodeXYZ np1( p1.node ), np2( p2.node );
865 RestrictLocalSize( ngMesh, 0.5*(np1+np2), (np1-np2).Modulus() );
868 cout << "Segment: " << seg.edgenr << " on SMESH face " << helper.GetMeshDS()->ShapeToIndex( face ) << endl
869 << "\tface index: " << seg.si << endl
870 << "\tp1: " << seg[0] << endl
871 << "\tp2: " << seg[1] << endl
872 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
873 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
874 //<< "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
875 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
876 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl;
877 //<< "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
881 if ( helper.GetPeriodicIndex() && 1 ) {
882 seg.epgeominfo[ 0 ].u = otherSeamParam;
883 seg.epgeominfo[ 1 ].u = otherSeamParam;
884 swap (seg.epgeominfo[0].v, seg.epgeominfo[1].v);
886 seg.epgeominfo[ 0 ].v = otherSeamParam;
887 seg.epgeominfo[ 1 ].v = otherSeamParam;
888 swap (seg.epgeominfo[0].u, seg.epgeominfo[1].u);
890 swap( seg[0], seg[1] );
891 swap( seg.epgeominfo[0].dist, seg.epgeominfo[1].dist );
892 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
893 ngMesh.AddSegment( seg );
895 cout << "Segment: " << seg.edgenr << endl
896 << "\t is SEAM (reverse) of the previous. "
897 << " Other " << (helper.GetPeriodicIndex() && 1 ? "U" : "V")
898 << " = " << otherSeamParam << endl;
901 else if ( fOri == TopAbs_INTERNAL )
903 swap( seg[0], seg[1] );
904 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
905 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
906 ngMesh.AddSegment( seg );
908 cout << "Segment: " << seg.edgenr << endl << "\t is REVERSE of the previous" << endl;
912 } // loop on geomEdge ancestors
914 if ( quadHelper ) // remember medium nodes of sub-meshes
916 SMDS_ElemIteratorPtr edges = smDS->GetElements();
917 while ( edges->more() )
919 const SMDS_MeshElement* e = edges->next();
920 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshEdge*>( e )))
926 } // case TopAbs_EDGE
928 case TopAbs_FACE: { // FACE
929 // ----------------------
930 const TopoDS_Face& geomFace = TopoDS::Face( sm->GetSubShape() );
931 helper.SetSubShape( geomFace );
932 bool isInternalFace = ( geomFace.Orientation() == TopAbs_INTERNAL );
934 // Find solids the geomFace bounds
935 int solidID1 = 0, solidID2 = 0;
936 StdMeshers_QuadToTriaAdaptor* quadAdaptor =
937 dynamic_cast<StdMeshers_QuadToTriaAdaptor*>( proxyMesh.get() );
940 solidID1 = occgeom.somap.FindIndex( quadAdaptor->GetShape() );
944 PShapeIteratorPtr solidIt = helper.GetAncestors( geomFace, *sm->GetFather(), TopAbs_SOLID);
945 while ( const TopoDS_Shape * solid = solidIt->next() )
947 int id = occgeom.somap.FindIndex ( *solid );
948 if ( solidID1 && id != solidID1 ) solidID2 = id;
952 // Add ng face descriptors of meshed faces
954 ngMesh.AddFaceDescriptor( netgen::FaceDescriptor( faceNgID, solidID1, solidID2, 0 ));
956 // if second oreder is required, even already meshed faces must be passed to NETGEN
957 int fID = occgeom.fmap.Add( geomFace );
958 if ( occgeom.facemeshstatus.Size() < fID ) occgeom.facemeshstatus.SetSize( fID );
959 occgeom.facemeshstatus[ fID-1 ] = netgen::FACE_MESHED_OK;
960 while ( fID < faceNgID ) // geomFace is already in occgeom.fmap, add a copy
962 fID = occgeom.fmap.Add( BRepBuilderAPI_Copy( geomFace, /*copyGeom=*/false ));
963 if ( occgeom.facemeshstatus.Size() < fID ) occgeom.facemeshstatus.SetSize( fID );
964 occgeom.facemeshstatus[ fID-1 ] = netgen::FACE_MESHED_OK;
966 // Problem with the second order in a quadrangular mesh remains.
967 // 1) All quadrangles generated by NETGEN are moved to an inexistent face
968 // by FillSMesh() (find "AddFaceDescriptor")
969 // 2) Temporary triangles generated by StdMeshers_QuadToTriaAdaptor
970 // are on faces where quadrangles were.
971 // Due to these 2 points, wrong geom faces are used while conversion to quadratic
972 // of the mentioned above quadrangles and triangles
974 // Orient the face correctly in solidID1 (issue 0020206)
975 bool reverse = false;
977 TopoDS_Shape solid = occgeom.somap( solidID1 );
978 TopAbs_Orientation faceOriInSolid = helper.GetSubShapeOri( solid, geomFace );
979 if ( faceOriInSolid >= 0 )
981 helper.IsReversedSubMesh( TopoDS::Face( geomFace.Oriented( faceOriInSolid )));
984 // Add surface elements
986 netgen::Element2d tri(3);
987 tri.SetIndex( faceNgID );
988 SMESH_TNodeXYZ xyz[3];
990 #ifdef DUMP_TRIANGLES
991 cout << "SMESH face " << helper.GetMeshDS()->ShapeToIndex( geomFace )
992 << " internal="<<isInternalFace << endl;
995 smDS = proxyMesh->GetSubMesh( geomFace );
997 SMDS_ElemIteratorPtr faces = smDS->GetElements();
998 while ( faces->more() )
1000 const SMDS_MeshElement* f = faces->next();
1001 if ( f->NbNodes() % 3 != 0 ) // not triangle
1003 PShapeIteratorPtr solidIt=helper.GetAncestors(geomFace,*sm->GetFather(),TopAbs_SOLID);
1004 if ( const TopoDS_Shape * solid = solidIt->next() )
1005 sm = _mesh->GetSubMesh( *solid );
1006 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1007 smError.reset( new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH,"Not triangle sub-mesh"));
1008 smError->myBadElements.push_back( f );
1012 for ( int i = 0; i < 3; ++i )
1014 const SMDS_MeshNode* node = f->GetNode( i ), * inFaceNode=0;
1017 // get node UV on face
1018 int shapeID = node->getshapeId();
1019 if ( helper.IsSeamShape( shapeID ))
1021 if ( helper.IsSeamShape( f->GetNodeWrap( i+1 )->getshapeId() ))
1022 inFaceNode = f->GetNodeWrap( i-1 );
1024 inFaceNode = f->GetNodeWrap( i+1 );
1026 gp_XY uv = helper.GetNodeUV( geomFace, node, inFaceNode );
1028 int ind = reverse ? 3-i : i+1;
1029 tri.GeomInfoPi(ind).u = uv.X();
1030 tri.GeomInfoPi(ind).v = uv.Y();
1031 tri.PNum (ind) = ngNodeId( node, ngMesh, nodeNgIdMap );
1034 // pass a triangle size to NG size-map
1035 double size = ( ( xyz[0] - xyz[1] ).Modulus() +
1036 ( xyz[1] - xyz[2] ).Modulus() +
1037 ( xyz[2] - xyz[0] ).Modulus() ) / 3;
1038 gp_XYZ gc = ( xyz[0] + xyz[1] + xyz[2] ) / 3;
1039 RestrictLocalSize( ngMesh, gc, size, /*overrideMinH=*/false );
1041 ngMesh.AddSurfaceElement (tri);
1042 #ifdef DUMP_TRIANGLES
1043 cout << tri << endl;
1046 if ( isInternalFace )
1048 swap( tri[1], tri[2] );
1049 ngMesh.AddSurfaceElement (tri);
1050 #ifdef DUMP_TRIANGLES
1051 cout << tri << endl;
1056 if ( quadHelper ) // remember medium nodes of sub-meshes
1058 SMDS_ElemIteratorPtr faces = smDS->GetElements();
1059 while ( faces->more() )
1061 const SMDS_MeshElement* f = faces->next();
1062 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshFace*>( f )))
1068 } // case TopAbs_FACE
1070 case TopAbs_VERTEX: { // VERTEX
1071 // --------------------------
1072 // issue 0021405. Add node only if a VERTEX is shared by a not meshed EDGE,
1073 // else netgen removes a free node and nodeVector becomes invalid
1074 PShapeIteratorPtr ansIt = helper.GetAncestors( sm->GetSubShape(),
1078 while ( const TopoDS_Shape* e = ansIt->next() )
1080 SMESH_subMesh* eSub = helper.GetMesh()->GetSubMesh( *e );
1081 if (( toAdd = ( eSub->IsEmpty() && !SMESH_Algo::isDegenerated( TopoDS::Edge( *e )))))
1086 SMDS_NodeIteratorPtr nodeIt = smDS->GetNodes();
1087 if ( nodeIt->more() )
1088 ngNodeId( nodeIt->next(), ngMesh, nodeNgIdMap );
1094 } // loop on submeshes
1097 nodeVec.resize( ngMesh.GetNP() + 1 );
1098 TNode2IdMap::iterator node_NgId, nodeNgIdEnd = nodeNgIdMap.end();
1099 for ( node_NgId = nodeNgIdMap.begin(); node_NgId != nodeNgIdEnd; ++node_NgId)
1100 nodeVec[ node_NgId->second ] = node_NgId->first;
1105 //================================================================================
1107 * \brief Duplicate mesh faces on internal geom faces
1109 //================================================================================
1111 void NETGENPlugin_Mesher::FixIntFaces(const netgen::OCCGeometry& occgeom,
1112 netgen::Mesh& ngMesh,
1113 NETGENPlugin_Internals& internalShapes)
1115 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1117 // find ng indices of internal faces
1119 for ( int ngFaceID = 1; ngFaceID <= occgeom.fmap.Extent(); ++ngFaceID )
1121 int smeshID = meshDS->ShapeToIndex( occgeom.fmap( ngFaceID ));
1122 if ( internalShapes.isInternalShape( smeshID ))
1123 ngFaceIds.insert( ngFaceID );
1125 if ( !ngFaceIds.empty() )
1128 int i, nbFaces = ngMesh.GetNSE();
1129 for ( i = 1; i <= nbFaces; ++i)
1131 netgen::Element2d elem = ngMesh.SurfaceElement(i);
1132 if ( ngFaceIds.count( elem.GetIndex() ))
1134 swap( elem[1], elem[2] );
1135 ngMesh.AddSurfaceElement (elem);
1141 //================================================================================
1143 * \brief Tries to heal the mesh on a FACE. The FACE is supposed to be partially
1144 * meshed due to NETGEN failure
1145 * \param [in] occgeom - geometry
1146 * \param [in,out] ngMesh - the mesh to fix
1147 * \param [inout] faceID - ID of the FACE to fix the mesh on
1148 * \return bool - is mesh is or becomes OK
1150 //================================================================================
1152 bool NETGENPlugin_Mesher::FixFaceMesh(const netgen::OCCGeometry& occgeom,
1153 netgen::Mesh& ngMesh,
1156 // we address a case where the FACE is almost fully meshed except small holes
1157 // of usually triangular shape at FACE boundary (IPAL52861)
1159 // The case appeared to be not simple: holes only look triangular but
1160 // indeed are a self intersecting polygon. A reason of the bug was in coincident
1161 // NG points on a seam edge. But the code below is very nice, leave it for
1166 if ( occgeom.fmap.Extent() < faceID )
1168 //const TopoDS_Face& face = TopoDS::Face( occgeom.fmap( faceID ));
1170 // find free links on the FACE
1171 NCollection_Map<Link> linkMap;
1172 for ( int iF = 1; iF <= ngMesh.GetNSE(); ++iF )
1174 const netgen::Element2d& elem = ngMesh.SurfaceElement(iF);
1175 if ( faceID != elem.GetIndex() )
1177 int n0 = elem[ elem.GetNP() - 1 ];
1178 for ( int i = 0; i < elem.GetNP(); ++i )
1181 Link link( n0, n1 );
1182 if ( !linkMap.Add( link ))
1183 linkMap.Remove( link );
1187 // add/remove boundary links
1188 for ( int iSeg = 1; iSeg <= ngMesh.GetNSeg(); ++iSeg )
1190 const netgen::Segment& seg = ngMesh.LineSegment( iSeg );
1191 if ( seg.si != faceID ) // !edgeIDs.Contains( seg.edgenr ))
1193 Link link( seg[1], seg[0] ); // reverse!!!
1194 if ( !linkMap.Add( link ))
1195 linkMap.Remove( link );
1197 if ( linkMap.IsEmpty() )
1199 if ( linkMap.Extent() < 3 )
1202 // make triangles of the links
1204 netgen::Element2d tri(3);
1205 tri.SetIndex ( faceID );
1207 NCollection_Map<Link>::Iterator linkIt( linkMap );
1208 Link link1 = linkIt.Value();
1209 // look for a link connected to link1
1210 NCollection_Map<Link>::Iterator linkIt2 = linkIt;
1211 for ( linkIt2.Next(); linkIt2.More(); linkIt2.Next() )
1213 const Link& link2 = linkIt2.Value();
1214 if ( link2.IsConnected( link1 ))
1216 // look for a link connected to both link1 and link2
1217 NCollection_Map<Link>::Iterator linkIt3 = linkIt2;
1218 for ( linkIt3.Next(); linkIt3.More(); linkIt3.Next() )
1220 const Link& link3 = linkIt3.Value();
1221 if ( link3.IsConnected( link1 ) &&
1222 link3.IsConnected( link2 ) )
1227 tri[2] = ( link2.Contains( link1.n1 ) ? link2.n1 : link3.n1 );
1228 if ( tri[0] == tri[2] || tri[1] == tri[2] )
1230 ngMesh.AddSurfaceElement( tri );
1232 // prepare for the next tria search
1233 if ( linkMap.Extent() == 3 )
1235 linkMap.Remove( link3 );
1236 linkMap.Remove( link2 );
1238 linkMap.Remove( link1 );
1239 link1 = linkIt.Value();
1252 //================================================================================
1253 // define gp_XY_Subtracted pointer to function calling gp_XY::Subtracted(gp_XY)
1254 gp_XY_FunPtr(Subtracted);
1255 //gp_XY_FunPtr(Added);
1257 //================================================================================
1259 * \brief Evaluate distance between two 2d points along the surface
1261 //================================================================================
1263 double evalDist( const gp_XY& uv1,
1265 const Handle(Geom_Surface)& surf,
1266 const int stopHandler=-1)
1268 if ( stopHandler > 0 ) // continue recursion
1270 gp_XY mid = SMESH_MesherHelper::GetMiddleUV( surf, uv1, uv2 );
1271 return evalDist( uv1,mid, surf, stopHandler-1 ) + evalDist( mid,uv2, surf, stopHandler-1 );
1273 double dist3D = surf->Value( uv1.X(), uv1.Y() ).Distance( surf->Value( uv2.X(), uv2.Y() ));
1274 if ( stopHandler == 0 ) // stop recursion
1277 // start recursion if necessary
1278 double dist2D = SMESH_MesherHelper::ApplyIn2D(surf, uv1, uv2, gp_XY_Subtracted, 0).Modulus();
1279 if ( fabs( dist3D - dist2D ) < dist2D * 1e-10 )
1280 return dist3D; // equal parametrization of a planar surface
1282 return evalDist( uv1, uv2, surf, 3 ); // start recursion
1285 //================================================================================
1287 * \brief Data of vertex internal in geom face
1289 //================================================================================
1293 gp_XY uv; //!< UV in face parametric space
1294 int ngId; //!< ng id of corrsponding node
1295 gp_XY uvClose; //!< UV of closest boundary node
1296 int ngIdClose; //!< ng id of closest boundary node
1299 //================================================================================
1301 * \brief Data of vertex internal in solid
1303 //================================================================================
1307 int ngId; //!< ng id of corresponding node
1308 int ngIdClose; //!< ng id of closest 2d mesh element
1309 int ngIdCloseN; //!< ng id of closest node of the closest 2d mesh element
1312 inline double dist2(const netgen::MeshPoint& p1, const netgen::MeshPoint& p2)
1314 return gp_Pnt( NGPOINT_COORDS(p1)).SquareDistance( gp_Pnt( NGPOINT_COORDS(p2)));
1318 //================================================================================
1320 * \brief Make netgen take internal vertices in faces into account by adding
1321 * segments including internal vertices
1323 * This function works in supposition that 1D mesh is already computed in ngMesh
1325 //================================================================================
1327 void NETGENPlugin_Mesher::AddIntVerticesInFaces(const netgen::OCCGeometry& occgeom,
1328 netgen::Mesh& ngMesh,
1329 vector<const SMDS_MeshNode*>& nodeVec,
1330 NETGENPlugin_Internals& internalShapes)
1332 if ((int) nodeVec.size() < ngMesh.GetNP() )
1333 nodeVec.resize( ngMesh.GetNP(), 0 );
1335 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1336 SMESH_MesherHelper helper( internalShapes.getMesh() );
1338 const map<int,list<int> >& face2Vert = internalShapes.getFacesWithVertices();
1339 map<int,list<int> >::const_iterator f2v = face2Vert.begin();
1340 for ( ; f2v != face2Vert.end(); ++f2v )
1342 const TopoDS_Face& face = TopoDS::Face( meshDS->IndexToShape( f2v->first ));
1343 if ( face.IsNull() ) continue;
1344 int faceNgID = occgeom.fmap.FindIndex (face);
1345 if ( faceNgID < 0 ) continue;
1347 TopLoc_Location loc;
1348 Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc);
1350 helper.SetSubShape( face );
1351 helper.SetElementsOnShape( true );
1353 // Get data of internal vertices and add them to ngMesh
1355 multimap< double, TIntVData > dist2VData; // sort vertices by distance from boundary nodes
1357 int i, nbSegInit = ngMesh.GetNSeg();
1359 // boundary characteristics
1360 double totSegLen2D = 0;
1363 const list<int>& iVertices = f2v->second;
1364 list<int>::const_iterator iv = iVertices.begin();
1365 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1368 // get node on vertex
1369 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1370 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1373 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1374 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1375 nV = SMESH_Algo::VertexNode( V, meshDS );
1376 if ( !nV ) continue;
1379 netgen::MeshPoint mp( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1380 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1381 vData.ngId = ngMesh.GetNP();
1382 nodeVec.push_back( nV );
1386 vData.uv = helper.GetNodeUV( face, nV, 0, &uvOK );
1387 if ( !uvOK ) helper.CheckNodeUV( face, nV, vData.uv, BRep_Tool::Tolerance(V),/*force=*/1);
1389 // loop on all segments of the face to find the node closest to vertex and to count
1390 // average segment 2d length
1391 double closeDist2 = numeric_limits<double>::max(), dist2;
1393 for (i = 1; i <= ngMesh.GetNSeg(); ++i)
1395 netgen::Segment & seg = ngMesh.LineSegment(i);
1396 if ( seg.si != faceNgID ) continue;
1398 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1400 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1401 if ( ngIdLast == seg[ iEnd ] ) continue;
1402 dist2 = helper.ApplyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1403 if ( dist2 < closeDist2 )
1404 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1405 ngIdLast = seg[ iEnd ];
1409 totSegLen2D += helper.ApplyIn2D(surf, uv[0], uv[1], gp_XY_Subtracted, false).Modulus();
1413 dist2VData.insert( make_pair( closeDist2, vData ));
1416 if ( totNbSeg == 0 ) break;
1417 double avgSegLen2d = totSegLen2D / totNbSeg;
1419 // Loop on vertices to add segments
1421 multimap< double, TIntVData >::iterator dist_vData = dist2VData.begin();
1422 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1424 double closeDist2 = dist_vData->first, dist2;
1425 TIntVData & vData = dist_vData->second;
1427 // try to find more close node among segments added for internal vertices
1428 for (i = nbSegInit+1; i <= ngMesh.GetNSeg(); ++i)
1430 netgen::Segment & seg = ngMesh.LineSegment(i);
1431 if ( seg.si != faceNgID ) continue;
1433 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1435 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1436 dist2 = helper.ApplyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1437 if ( dist2 < closeDist2 )
1438 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1441 // decide whether to use the closest node as the second end of segment or to
1442 // create a new point
1443 int segEnd1 = vData.ngId;
1444 int segEnd2 = vData.ngIdClose; // to use closest node
1445 gp_XY uvV = vData.uv, uvP = vData.uvClose;
1446 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1447 double nodeDist2D = sqrt( closeDist2 );
1448 double nodeDist3D = evalDist( vData.uv, vData.uvClose, surf );
1449 bool avgLenOK = ( avgSegLen2d < 0.75 * nodeDist2D );
1450 bool hintLenOK = ( segLenHint < 0.75 * nodeDist3D );
1451 //cout << "uvV " << uvV.X() <<","<<uvV.Y() << " ";
1452 if ( hintLenOK || avgLenOK )
1454 // create a point between the closest node and V
1457 double r = min( 0.5, ( hintLenOK ? segLenHint/nodeDist3D : avgSegLen2d/nodeDist2D ));
1458 // direction from V to closet node in 2D
1459 gp_Dir2d v2n( helper.ApplyIn2D(surf, uvP, uvV, gp_XY_Subtracted, false ));
1461 uvP = vData.uv + r * nodeDist2D * v2n.XY();
1462 gp_Pnt P = surf->Value( uvP.X(), uvP.Y() ).Transformed( loc );
1464 netgen::MeshPoint mp( netgen::Point<3> (P.X(), P.Y(), P.Z()));
1465 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1466 segEnd2 = ngMesh.GetNP();
1467 //cout << "Middle " << r << " uv " << uvP.X() << "," << uvP.Y() << "( " << ngMesh.Point(segEnd2).X()<<","<<ngMesh.Point(segEnd2).Y()<<","<<ngMesh.Point(segEnd2).Z()<<" )"<< endl;
1468 SMDS_MeshNode * nP = helper.AddNode(P.X(), P.Y(), P.Z());
1469 nodeVec.push_back( nP );
1471 //else cout << "at Node " << " uv " << uvP.X() << "," << uvP.Y() << endl;
1474 netgen::Segment seg;
1476 if ( segEnd1 > segEnd2 ) swap( segEnd1, segEnd2 ), swap( uvV, uvP );
1477 seg[0] = segEnd1; // ng node id
1478 seg[1] = segEnd2; // ng node id
1479 seg.edgenr = ngMesh.GetNSeg() + 1;// segment id
1482 seg.epgeominfo[ 0 ].dist = 0; // param on curve
1483 seg.epgeominfo[ 0 ].u = uvV.X();
1484 seg.epgeominfo[ 0 ].v = uvV.Y();
1485 seg.epgeominfo[ 1 ].dist = 1; // param on curve
1486 seg.epgeominfo[ 1 ].u = uvP.X();
1487 seg.epgeominfo[ 1 ].v = uvP.Y();
1489 // seg.epgeominfo[ 0 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1490 // seg.epgeominfo[ 1 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1492 ngMesh.AddSegment (seg);
1494 // add reverse segment
1495 swap( seg[0], seg[1] );
1496 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1497 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1498 ngMesh.AddSegment (seg);
1504 //================================================================================
1506 * \brief Make netgen take internal vertices in solids into account by adding
1507 * faces including internal vertices
1509 * This function works in supposition that 2D mesh is already computed in ngMesh
1511 //================================================================================
1513 void NETGENPlugin_Mesher::AddIntVerticesInSolids(const netgen::OCCGeometry& occgeom,
1514 netgen::Mesh& ngMesh,
1515 vector<const SMDS_MeshNode*>& nodeVec,
1516 NETGENPlugin_Internals& internalShapes)
1518 #ifdef DUMP_TRIANGLES_SCRIPT
1519 // create a python script making a mesh containing triangles added for internal vertices
1520 ofstream py(DUMP_TRIANGLES_SCRIPT);
1521 py << "import SMESH"<< endl
1522 << "from salome.smesh import smeshBuilder"<<endl
1523 << "smesh = smeshBuilder.New(salome.myStudy)"<<endl
1524 << "m = smesh.Mesh(name='triangles')" << endl;
1526 if ((int) nodeVec.size() < ngMesh.GetNP() )
1527 nodeVec.resize( ngMesh.GetNP(), 0 );
1529 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1530 SMESH_MesherHelper helper( internalShapes.getMesh() );
1532 const map<int,list<int> >& so2Vert = internalShapes.getSolidsWithVertices();
1533 map<int,list<int> >::const_iterator s2v = so2Vert.begin();
1534 for ( ; s2v != so2Vert.end(); ++s2v )
1536 const TopoDS_Shape& solid = meshDS->IndexToShape( s2v->first );
1537 if ( solid.IsNull() ) continue;
1538 int solidNgID = occgeom.somap.FindIndex (solid);
1539 if ( solidNgID < 0 && !occgeom.somap.IsEmpty() ) continue;
1541 helper.SetSubShape( solid );
1542 helper.SetElementsOnShape( true );
1544 // find ng indices of faces within the solid
1546 for (TopExp_Explorer fExp(solid, TopAbs_FACE); fExp.More(); fExp.Next() )
1547 ngFaceIds.insert( occgeom.fmap.FindIndex( fExp.Current() ));
1548 if ( ngFaceIds.size() == 1 && *ngFaceIds.begin() == 0 )
1549 ngFaceIds.insert( 1 );
1551 // Get data of internal vertices and add them to ngMesh
1553 multimap< double, TIntVSoData > dist2VData; // sort vertices by distance from ng faces
1555 int i, nbFaceInit = ngMesh.GetNSE();
1557 // boundary characteristics
1558 double totSegLen = 0;
1561 const list<int>& iVertices = s2v->second;
1562 list<int>::const_iterator iv = iVertices.begin();
1563 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1566 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1568 // get node on vertex
1569 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1572 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1573 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1574 nV = SMESH_Algo::VertexNode( V, meshDS );
1575 if ( !nV ) continue;
1578 netgen::MeshPoint mpV( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1579 ngMesh.AddPoint ( mpV, 1, netgen::FIXEDPOINT );
1580 vData.ngId = ngMesh.GetNP();
1581 nodeVec.push_back( nV );
1583 // loop on all 2d elements to find the one closest to vertex and to count
1584 // average segment length
1585 double closeDist2 = numeric_limits<double>::max(), avgDist2;
1586 for (i = 1; i <= ngMesh.GetNSE(); ++i)
1588 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1589 if ( !ngFaceIds.count( elem.GetIndex() )) continue;
1591 multimap< double, int> dist2nID; // sort nodes of element by distance from V
1592 for ( int j = 0; j < elem.GetNP(); ++j)
1594 netgen::MeshPoint mp = ngMesh.Point( elem[j] );
1595 double d2 = dist2( mpV, mp );
1596 dist2nID.insert( make_pair( d2, elem[j] ));
1597 avgDist2 += d2 / elem.GetNP();
1599 totNbSeg++, totSegLen+= sqrt( dist2( mp, ngMesh.Point( elem[(j+1)%elem.GetNP()])));
1601 double dist = dist2nID.begin()->first; //avgDist2;
1602 if ( dist < closeDist2 )
1603 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= dist;
1605 dist2VData.insert( make_pair( closeDist2, vData ));
1608 if ( totNbSeg == 0 ) break;
1609 double avgSegLen = totSegLen / totNbSeg;
1611 // Loop on vertices to add triangles
1613 multimap< double, TIntVSoData >::iterator dist_vData = dist2VData.begin();
1614 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1616 double closeDist2 = dist_vData->first;
1617 TIntVSoData & vData = dist_vData->second;
1619 const netgen::MeshPoint& mpV = ngMesh.Point( vData.ngId );
1621 // try to find more close face among ones added for internal vertices
1622 for (i = nbFaceInit+1; i <= ngMesh.GetNSE(); ++i)
1624 double avgDist2 = 0;
1625 multimap< double, int> dist2nID;
1626 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1627 for ( int j = 0; j < elem.GetNP(); ++j)
1629 double d = dist2( mpV, ngMesh.Point( elem[j] ));
1630 dist2nID.insert( make_pair( d, elem[j] ));
1631 avgDist2 += d / elem.GetNP();
1632 if ( avgDist2 < closeDist2 )
1633 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= avgDist2;
1636 // sort nodes of the closest face by angle with vector from V to the closest node
1637 const double tol = numeric_limits<double>::min();
1638 map< double, int > angle2ID;
1639 const netgen::Element2d& closeFace = ngMesh.SurfaceElement( vData.ngIdClose );
1640 netgen::MeshPoint mp[2];
1641 mp[0] = ngMesh.Point( vData.ngIdCloseN );
1642 gp_XYZ p1( NGPOINT_COORDS( mp[0] ));
1643 gp_XYZ pV( NGPOINT_COORDS( mpV ));
1644 gp_Vec v2p1( pV, p1 );
1645 double distN1 = v2p1.Magnitude();
1646 if ( distN1 <= tol ) continue;
1648 for ( int j = 0; j < closeFace.GetNP(); ++j)
1650 mp[1] = ngMesh.Point( closeFace[j] );
1651 gp_Vec v2p( pV, gp_Pnt( NGPOINT_COORDS( mp[1] )) );
1652 angle2ID.insert( make_pair( v2p1.Angle( v2p ), closeFace[j]));
1654 // get node with angle of 60 degrees or greater
1655 map< double, int >::iterator angle_id = angle2ID.lower_bound( 60. * M_PI / 180. );
1656 if ( angle_id == angle2ID.end() ) angle_id = --angle2ID.end();
1657 const double minAngle = 30. * M_PI / 180.;
1658 const double angle = angle_id->first;
1659 bool angleOK = ( angle > minAngle );
1661 // find points to create a triangle
1662 netgen::Element2d tri(3);
1664 tri[0] = vData.ngId;
1665 tri[1] = vData.ngIdCloseN; // to use the closest nodes
1666 tri[2] = angle_id->second; // to use the node with best angle
1668 // decide whether to use the closest node and the node with best angle or to create new ones
1669 for ( int isBestAngleN = 0; isBestAngleN < 2; ++isBestAngleN )
1671 bool createNew = !angleOK; //, distOK = true;
1673 int triInd = isBestAngleN ? 2 : 1;
1674 mp[isBestAngleN] = ngMesh.Point( tri[triInd] );
1679 double distN2 = sqrt( dist2( mpV, mp[isBestAngleN]));
1680 createNew = ( fabs( distN2 - distN1 ) > 0.25 * distN1 );
1682 else if ( angle < tol )
1684 v2p1.SetX( v2p1.X() + 1e-3 );
1690 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1691 bool avgLenOK = ( avgSegLen < 0.75 * distN1 );
1692 bool hintLenOK = ( segLenHint < 0.75 * distN1 );
1693 createNew = (createNew || avgLenOK || hintLenOK );
1694 // we create a new node not closer than 0.5 to the closest face
1695 // in order not to clash with other close face
1696 double r = min( 0.5, ( hintLenOK ? segLenHint : avgSegLen ) / distN1 );
1697 distFromV = r * distN1;
1701 // create a new point, between the node and the vertex if angleOK
1702 gp_XYZ p( NGPOINT_COORDS( mp[isBestAngleN] ));
1703 gp_Vec v2p( pV, p ); v2p.Normalize();
1704 if ( isBestAngleN && !angleOK )
1705 p = p1 + gp_Dir( v2p.XYZ() - v2p1.XYZ()).XYZ() * distN1 * 0.95;
1707 p = pV + v2p.XYZ() * distFromV;
1709 if ( !isBestAngleN ) p1 = p, distN1 = distFromV;
1711 mp[isBestAngleN].SetPoint( netgen::Point<3> (p.X(), p.Y(), p.Z()));
1712 ngMesh.AddPoint ( mp[isBestAngleN], 1, netgen::SURFACEPOINT );
1713 tri[triInd] = ngMesh.GetNP();
1714 nodeVec.push_back( helper.AddNode( p.X(), p.Y(), p.Z()) );
1717 ngMesh.AddSurfaceElement (tri);
1718 swap( tri[1], tri[2] );
1719 ngMesh.AddSurfaceElement (tri);
1721 #ifdef DUMP_TRIANGLES_SCRIPT
1722 py << "n1 = m.AddNode( "<< mpV(0)<<", "<< mpV(1)<<", "<< mpV(2)<<") "<< endl
1723 << "n2 = m.AddNode( "<< mp[0](0)<<", "<< mp[0](1)<<", "<< mp[0](2)<<") "<< endl
1724 << "n3 = m.AddNode( "<< mp[1](0)<<", "<< mp[1](1)<<", "<< mp[1](2)<<" )" << endl
1725 << "m.AddFace([n1,n2,n3])" << endl;
1727 } // loop on internal vertices of a solid
1729 } // loop on solids with internal vertices
1732 //================================================================================
1734 * \brief Fill netgen mesh with segments of a FACE
1735 * \param ngMesh - netgen mesh
1736 * \param geom - container of OCCT geometry to mesh
1737 * \param wires - data of nodes on FACE boundary
1738 * \param helper - mesher helper holding the FACE
1739 * \param nodeVec - vector of nodes in which node index == netgen ID
1740 * \retval SMESH_ComputeErrorPtr - error description
1742 //================================================================================
1744 SMESH_ComputeErrorPtr
1745 NETGENPlugin_Mesher::AddSegmentsToMesh(netgen::Mesh& ngMesh,
1746 netgen::OCCGeometry& geom,
1747 const TSideVector& wires,
1748 SMESH_MesherHelper& helper,
1749 vector< const SMDS_MeshNode* > & nodeVec,
1750 const bool overrideMinH)
1752 // ----------------------------
1753 // Check wires and count nodes
1754 // ----------------------------
1756 for ( size_t iW = 0; iW < wires.size(); ++iW )
1758 StdMeshers_FaceSidePtr wire = wires[ iW ];
1759 if ( wire->MissVertexNode() )
1761 // Commented for issue 0020960. It worked for the case, let's wait for case where it doesn't.
1762 // It seems that there is no reason for this limitation
1764 // (new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH, "Missing nodes on vertices"));
1766 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1767 if ((int) uvPtVec.size() != wire->NbPoints() )
1768 return SMESH_ComputeError::New(COMPERR_BAD_INPUT_MESH,
1769 SMESH_Comment("Unexpected nb of points on wire ") << iW
1770 << ": " << uvPtVec.size()<<" != "<<wire->NbPoints());
1771 nbNodes += wire->NbPoints();
1773 nodeVec.reserve( nodeVec.size() + nbNodes + 1 );
1774 if ( nodeVec.empty() )
1775 nodeVec.push_back( 0 );
1777 // -----------------
1779 // -----------------
1781 const bool wasNgMeshEmpty = ( ngMesh.GetNP() < 1 ); /* true => this method is called by
1782 NETGENPlugin_NETGEN_2D_ONLY */
1784 // map for nodes on vertices since they can be shared between wires
1785 // ( issue 0020676, face_int_box.brep) and nodes built by NETGEN
1786 map<const SMDS_MeshNode*, int > node2ngID;
1787 if ( !wasNgMeshEmpty ) // fill node2ngID with nodes built by NETGEN
1789 set< int > subIDs; // ids of sub-shapes of the FACE
1790 for ( size_t iW = 0; iW < wires.size(); ++iW )
1792 StdMeshers_FaceSidePtr wire = wires[ iW ];
1793 for ( int iE = 0, nbE = wire->NbEdges(); iE < nbE; ++iE )
1795 subIDs.insert( wire->EdgeID( iE ));
1796 subIDs.insert( helper.GetMeshDS()->ShapeToIndex( wire->FirstVertex( iE )));
1799 for ( size_t ngID = 1; ngID < nodeVec.size(); ++ngID )
1800 if ( subIDs.count( nodeVec[ngID]->getshapeId() ))
1801 node2ngID.insert( make_pair( nodeVec[ngID], ngID ));
1804 const int solidID = 0, faceID = geom.fmap.FindIndex( helper.GetSubShape() );
1805 if ( ngMesh.GetNFD() < 1 )
1806 ngMesh.AddFaceDescriptor( netgen::FaceDescriptor( faceID, solidID, solidID, 0 ));
1808 for ( size_t iW = 0; iW < wires.size(); ++iW )
1810 StdMeshers_FaceSidePtr wire = wires[ iW ];
1811 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1812 const int nbSegments = wire->NbPoints() - 1;
1814 // assure the 1st node to be in node2ngID, which is needed to correctly
1815 // "close chain of segments" (see below) in case if the 1st node is not
1816 // onVertex because it is on a Viscous layer
1817 node2ngID.insert( make_pair( uvPtVec[ 0 ].node, ngMesh.GetNP() + 1 ));
1819 // compute length of every segment
1820 vector<double> segLen( nbSegments );
1821 for ( int i = 0; i < nbSegments; ++i )
1822 segLen[i] = SMESH_TNodeXYZ( uvPtVec[ i ].node ).Distance( uvPtVec[ i+1 ].node );
1824 int edgeID = 1, posID = -2;
1825 bool isInternalWire = false;
1826 double vertexNormPar = 0;
1827 //const int prevNbNGSeg = ngMesh.GetNSeg();
1828 for ( int i = 0; i < nbSegments; ++i ) // loop on segments
1830 // Add the first point of a segment
1832 const SMDS_MeshNode * n = uvPtVec[ i ].node;
1833 const int posShapeID = n->getshapeId();
1834 bool onVertex = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX );
1835 bool onEdge = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_EDGE );
1837 // skip nodes on degenerated edges
1838 if ( helper.IsDegenShape( posShapeID ) &&
1839 helper.IsDegenShape( uvPtVec[ i+1 ].node->getshapeId() ))
1842 int ngID1 = ngMesh.GetNP() + 1, ngID2 = ngID1+1;
1843 if ( onVertex || ( !wasNgMeshEmpty && onEdge ) || helper.IsRealSeam( posShapeID ))
1844 ngID1 = node2ngID.insert( make_pair( n, ngID1 )).first->second;
1845 if ( ngID1 > ngMesh.GetNP() )
1847 netgen::MeshPoint mp( netgen::Point<3> (n->X(), n->Y(), n->Z()) );
1848 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1849 nodeVec.push_back( n );
1851 else // n is in ngMesh already, and ngID2 in prev segment is wrong
1853 ngID2 = ngMesh.GetNP() + 1;
1854 if ( i > 0 ) // prev segment belongs to same wire
1856 netgen::Segment& prevSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
1863 netgen::Segment seg;
1865 seg[0] = ngID1; // ng node id
1866 seg[1] = ngID2; // ng node id
1867 seg.edgenr = ngMesh.GetNSeg() + 1; // ng segment id
1868 seg.si = faceID; // = geom.fmap.FindIndex (face);
1870 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1872 const UVPtStruct& pnt = uvPtVec[ i + iEnd ];
1874 seg.epgeominfo[ iEnd ].dist = pnt.param; // param on curve
1875 seg.epgeominfo[ iEnd ].u = pnt.u;
1876 seg.epgeominfo[ iEnd ].v = pnt.v;
1878 // find out edge id and node parameter on edge
1879 onVertex = ( pnt.normParam + 1e-10 > vertexNormPar );
1880 if ( onVertex || posShapeID != posID )
1883 double normParam = pnt.normParam;
1885 normParam = 0.5 * ( uvPtVec[ i ].normParam + uvPtVec[ i+1 ].normParam );
1886 int edgeIndexInWire = wire->EdgeIndex( normParam );
1887 vertexNormPar = wire->LastParameter( edgeIndexInWire );
1888 const TopoDS_Edge& edge = wire->Edge( edgeIndexInWire );
1889 edgeID = geom.emap.FindIndex( edge );
1891 isInternalWire = ( edge.Orientation() == TopAbs_INTERNAL );
1892 // if ( onVertex ) // param on curve is different on each of two edges
1893 // seg.epgeominfo[ iEnd ].dist = helper.GetNodeU( edge, pnt.node );
1895 seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
1898 ngMesh.AddSegment (seg);
1900 // restrict size of elements near the segment
1901 SMESH_TNodeXYZ np1( n ), np2( uvPtVec[ i+1 ].node );
1902 // get an average size of adjacent segments to avoid sharp change of
1903 // element size (regression on issue 0020452, note 0010898)
1904 int iPrev = SMESH_MesherHelper::WrapIndex( i-1, nbSegments );
1905 int iNext = SMESH_MesherHelper::WrapIndex( i+1, nbSegments );
1906 double sumH = segLen[ iPrev ] + segLen[ i ] + segLen[ iNext ];
1907 int nbSeg = ( int( segLen[ iPrev ] > sumH / 100.) +
1908 int( segLen[ i ] > sumH / 100.) +
1909 int( segLen[ iNext ] > sumH / 100.));
1911 RestrictLocalSize( ngMesh, 0.5*(np1+np2), sumH / nbSeg, overrideMinH );
1913 if ( isInternalWire )
1915 swap (seg[0], seg[1]);
1916 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1917 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1918 ngMesh.AddSegment (seg);
1920 } // loop on segments on a wire
1922 // close chain of segments
1923 if ( nbSegments > 0 )
1925 netgen::Segment& lastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() - int( isInternalWire ));
1926 const SMDS_MeshNode * lastNode = uvPtVec.back().node;
1927 lastSeg[1] = node2ngID.insert( make_pair( lastNode, lastSeg[1] )).first->second;
1928 if ( lastSeg[1] > ngMesh.GetNP() )
1930 netgen::MeshPoint mp( netgen::Point<3> (lastNode->X(), lastNode->Y(), lastNode->Z()) );
1931 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1932 nodeVec.push_back( lastNode );
1934 if ( isInternalWire )
1936 netgen::Segment& realLastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
1937 realLastSeg[0] = lastSeg[1];
1941 #ifdef DUMP_SEGMENTS
1942 cout << "BEGIN WIRE " << iW << endl;
1943 for ( int i = prevNbNGSeg+1; i <= ngMesh.GetNSeg(); ++i )
1945 netgen::Segment& seg = ngMesh.LineSegment( i );
1947 netgen::Segment& prevSeg = ngMesh.LineSegment( i-1 );
1948 if ( seg[0] == prevSeg[1] && seg[1] == prevSeg[0] )
1950 cout << "Segment: " << seg.edgenr << endl << "\tis REVRESE of the previous one" << endl;
1954 cout << "Segment: " << seg.edgenr << endl
1955 << "\tp1: " << seg[0] << " n" << nodeVec[ seg[0]]->GetID() << endl
1956 << "\tp2: " << seg[1] << " n" << nodeVec[ seg[1]]->GetID() << endl
1957 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
1958 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
1959 << "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
1960 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
1961 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl
1962 << "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
1964 cout << "--END WIRE " << iW << endl;
1967 } // loop on WIREs of a FACE
1969 // add a segment instead of an internal vertex
1970 if ( wasNgMeshEmpty )
1972 NETGENPlugin_Internals intShapes( *helper.GetMesh(), helper.GetSubShape(), /*is3D=*/false );
1973 AddIntVerticesInFaces( geom, ngMesh, nodeVec, intShapes );
1975 ngMesh.CalcSurfacesOfNode();
1980 //================================================================================
1982 * \brief Fill SMESH mesh according to contents of netgen mesh
1983 * \param occgeo - container of OCCT geometry to mesh
1984 * \param ngMesh - netgen mesh
1985 * \param initState - bn of entities in netgen mesh before computing
1986 * \param sMesh - SMESH mesh to fill in
1987 * \param nodeVec - vector of nodes in which node index == netgen ID
1988 * \param comment - returns problem description
1989 * \param quadHelper - holder of medium nodes of sub-meshes
1990 * \retval int - error
1992 //================================================================================
1994 int NETGENPlugin_Mesher::FillSMesh(const netgen::OCCGeometry& occgeo,
1995 netgen::Mesh& ngMesh,
1996 const NETGENPlugin_ngMeshInfo& initState,
1998 std::vector<const SMDS_MeshNode*>& nodeVec,
1999 SMESH_Comment& comment,
2000 SMESH_MesherHelper* quadHelper)
2002 int nbNod = ngMesh.GetNP();
2003 int nbSeg = ngMesh.GetNSeg();
2004 int nbFac = ngMesh.GetNSE();
2005 int nbVol = ngMesh.GetNE();
2007 SMESHDS_Mesh* meshDS = sMesh.GetMeshDS();
2009 // quadHelper is used for either
2010 // 1) making quadratic elements when a lower dimention mesh is loaded
2011 // to SMESH before convertion to quadratic by NETGEN
2012 // 2) sewing of quadratic elements with quadratic elements of sub-meshes
2013 if ( quadHelper && !quadHelper->GetIsQuadratic() && quadHelper->GetTLinkNodeMap().empty() )
2016 // -------------------------------------
2017 // Create and insert nodes into nodeVec
2018 // -------------------------------------
2020 nodeVec.resize( nbNod + 1 );
2021 int i, nbInitNod = initState._nbNodes;
2022 for (i = nbInitNod+1; i <= nbNod; ++i )
2024 const netgen::MeshPoint& ngPoint = ngMesh.Point(i);
2025 SMDS_MeshNode* node = NULL;
2026 TopoDS_Vertex aVert;
2027 // First, netgen creates nodes on vertices in occgeo.vmap,
2028 // so node index corresponds to vertex index
2029 // but (issue 0020776) netgen does not create nodes with equal coordinates
2030 if ( i-nbInitNod <= occgeo.vmap.Extent() )
2032 gp_Pnt p ( NGPOINT_COORDS(ngPoint) );
2033 for (int iV = i-nbInitNod; aVert.IsNull() && iV <= occgeo.vmap.Extent(); ++iV)
2035 aVert = TopoDS::Vertex( occgeo.vmap( iV ));
2036 gp_Pnt pV = BRep_Tool::Pnt( aVert );
2037 if ( p.SquareDistance( pV ) > 1e-20 )
2040 node = const_cast<SMDS_MeshNode*>( SMESH_Algo::VertexNode( aVert, meshDS ));
2043 if (!node) // node not found on vertex
2045 node = meshDS->AddNode( NGPOINT_COORDS( ngPoint ));
2046 if (!aVert.IsNull())
2047 meshDS->SetNodeOnVertex(node, aVert);
2052 // -------------------------------------------
2053 // Create mesh segments along geometric edges
2054 // -------------------------------------------
2056 int nbInitSeg = initState._nbSegments;
2057 for (i = nbInitSeg+1; i <= nbSeg; ++i )
2059 const netgen::Segment& seg = ngMesh.LineSegment(i);
2061 int pinds[3] = { seg.pnums[0], seg.pnums[1], seg.pnums[2] };
2064 for (int j=0; j < 3; ++j)
2066 int pind = pinds[j];
2067 if (pind <= 0 || !nodeVec_ACCESS(pind))
2075 int aGeomEdgeInd = seg.epgeominfo[j].edgenr;
2076 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
2077 aEdge = TopoDS::Edge(occgeo.emap(aGeomEdgeInd));
2079 param = seg.epgeominfo[j].dist;
2082 else // middle point
2084 param = param2 * 0.5;
2086 if (!aEdge.IsNull() && nodeVec_ACCESS(pind)->getshapeId() < 1)
2088 meshDS->SetNodeOnEdge(nodeVec_ACCESS(pind), aEdge, param);
2093 SMDS_MeshEdge* edge = 0;
2094 if (nbp == 2) // second order ?
2096 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1])))
2098 if ( quadHelper ) // final mesh must be quadratic
2099 edge = quadHelper->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
2101 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
2105 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
2106 nodeVec_ACCESS(pinds[2])))
2108 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
2109 nodeVec_ACCESS(pinds[2]));
2113 if ( comment.empty() ) comment << "Cannot create a mesh edge";
2114 MESSAGE("Cannot create a mesh edge");
2115 nbSeg = nbFac = nbVol = 0;
2118 if ( !aEdge.IsNull() && edge->getshapeId() < 1 )
2119 meshDS->SetMeshElementOnShape(edge, aEdge);
2121 else if ( comment.empty() )
2123 comment << "Invalid netgen segment #" << i;
2127 // ----------------------------------------
2128 // Create mesh faces along geometric faces
2129 // ----------------------------------------
2131 int nbInitFac = initState._nbFaces;
2132 int quadFaceID = ngMesh.GetNFD() + 1;
2133 if ( nbInitFac < nbFac )
2134 // add a faces descriptor to exclude qudrangle elements generated by NETGEN
2135 // from computation of 3D mesh
2136 ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(quadFaceID, /*solid1=*/0, /*solid2=*/0, 0));
2138 vector<const SMDS_MeshNode*> nodes;
2139 for (i = nbInitFac+1; i <= nbFac; ++i )
2141 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
2142 const int aGeomFaceInd = elem.GetIndex();
2144 if (aGeomFaceInd > 0 && aGeomFaceInd <= occgeo.fmap.Extent())
2145 aFace = TopoDS::Face(occgeo.fmap(aGeomFaceInd));
2147 for ( int j = 1; j <= elem.GetNP(); ++j )
2149 int pind = elem.PNum(j);
2150 if ( pind < 1 || pind >= (int) nodeVec.size() )
2152 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind))
2154 nodes.push_back( node );
2155 if (!aFace.IsNull() && node->getshapeId() < 1)
2157 const netgen::PointGeomInfo& pgi = elem.GeomInfoPi(j);
2158 meshDS->SetNodeOnFace(node, aFace, pgi.u, pgi.v);
2162 if ((int) nodes.size() != elem.GetNP() )
2164 if ( comment.empty() )
2165 comment << "Invalid netgen 2d element #" << i;
2166 continue; // bad node ids
2168 SMDS_MeshFace* face = NULL;
2169 switch (elem.GetType())
2172 if ( quadHelper ) // final mesh must be quadratic
2173 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2]);
2175 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2]);
2178 if ( quadHelper ) // final mesh must be quadratic
2179 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2181 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2182 // exclude qudrangle elements from computation of 3D mesh
2183 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2186 nodes[5] = mediumNode( nodes[0],nodes[1],nodes[5], quadHelper );
2187 nodes[3] = mediumNode( nodes[1],nodes[2],nodes[3], quadHelper );
2188 nodes[4] = mediumNode( nodes[2],nodes[0],nodes[4], quadHelper );
2189 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[5],nodes[3],nodes[4]);
2192 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2193 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2194 nodes[5] = mediumNode( nodes[2],nodes[3],nodes[5], quadHelper );
2195 nodes[6] = mediumNode( nodes[3],nodes[0],nodes[6], quadHelper );
2196 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3],
2197 nodes[4],nodes[7],nodes[5],nodes[6]);
2198 // exclude qudrangle elements from computation of 3D mesh
2199 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2202 MESSAGE("NETGEN created a face of unexpected type, ignoring");
2207 if ( comment.empty() ) comment << "Cannot create a mesh face";
2208 MESSAGE("Cannot create a mesh face");
2209 nbSeg = nbFac = nbVol = 0;
2212 if ( !aFace.IsNull() )
2213 meshDS->SetMeshElementOnShape( face, aFace );
2216 // ------------------
2217 // Create tetrahedra
2218 // ------------------
2220 for ( i = 1; i <= nbVol; ++i )
2222 const netgen::Element& elem = ngMesh.VolumeElement(i);
2223 int aSolidInd = elem.GetIndex();
2224 TopoDS_Solid aSolid;
2225 if ( aSolidInd > 0 && aSolidInd <= occgeo.somap.Extent() )
2226 aSolid = TopoDS::Solid(occgeo.somap(aSolidInd));
2228 for ( int j = 1; j <= elem.GetNP(); ++j )
2230 int pind = elem.PNum(j);
2231 if ( pind < 1 || pind >= (int)nodeVec.size() )
2233 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind) )
2235 nodes.push_back(node);
2236 if ( !aSolid.IsNull() && node->getshapeId() < 1 )
2237 meshDS->SetNodeInVolume(node, aSolid);
2240 if ((int) nodes.size() != elem.GetNP() )
2242 if ( comment.empty() )
2243 comment << "Invalid netgen 3d element #" << i;
2246 SMDS_MeshVolume* vol = NULL;
2247 switch ( elem.GetType() )
2250 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3]);
2253 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2254 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2255 nodes[5] = mediumNode( nodes[2],nodes[0],nodes[5], quadHelper );
2256 nodes[6] = mediumNode( nodes[0],nodes[3],nodes[6], quadHelper );
2257 nodes[8] = mediumNode( nodes[1],nodes[3],nodes[8], quadHelper );
2258 nodes[9] = mediumNode( nodes[2],nodes[3],nodes[9], quadHelper );
2259 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3],
2260 nodes[4],nodes[7],nodes[5],nodes[6],nodes[8],nodes[9]);
2263 MESSAGE("NETGEN created a volume of unexpected type, ignoring");
2268 if ( comment.empty() ) comment << "Cannot create a mesh volume";
2269 MESSAGE("Cannot create a mesh volume");
2270 nbSeg = nbFac = nbVol = 0;
2273 if (!aSolid.IsNull())
2274 meshDS->SetMeshElementOnShape(vol, aSolid);
2276 return comment.empty() ? 0 : 1;
2281 //================================================================================
2283 * \brief Restrict size of elements on the given edge
2285 //================================================================================
2287 void setLocalSize(const TopoDS_Edge& edge,
2291 if ( size <= std::numeric_limits<double>::min() )
2293 Standard_Real u1, u2;
2294 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, u1, u2);
2295 if ( curve.IsNull() )
2297 TopoDS_Iterator vIt( edge );
2298 if ( !vIt.More() ) return;
2299 gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex( vIt.Value() ));
2300 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
2304 const int nb = (int)( 1.5 * SMESH_Algo::EdgeLength( edge ) / size );
2305 Standard_Real delta = (u2-u1)/nb;
2306 for(int i=0; i<nb; i++)
2308 Standard_Real u = u1 + delta*i;
2309 gp_Pnt p = curve->Value(u);
2310 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
2311 netgen::Point3d pi(p.X(), p.Y(), p.Z());
2312 double resultSize = mesh.GetH(pi);
2313 if ( resultSize - size > 0.1*size )
2314 // netgen does restriction iff oldH/newH > 1.2 (localh.cpp:136)
2315 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), resultSize/1.201 );
2320 //================================================================================
2322 * \brief Convert error into text
2324 //================================================================================
2326 std::string text(int err)
2331 SMESH_Comment("Error in netgen::OCCGenerateMesh() at ") << netgen::multithread.task;
2334 //================================================================================
2336 * \brief Convert exception into text
2338 //================================================================================
2340 std::string text(Standard_Failure& ex)
2342 SMESH_Comment str("Exception in netgen::OCCGenerateMesh()");
2343 str << " at " << netgen::multithread.task
2344 << ": " << ex.DynamicType()->Name();
2345 if ( ex.GetMessageString() && strlen( ex.GetMessageString() ))
2346 str << ": " << ex.GetMessageString();
2349 //================================================================================
2351 * \brief Convert exception into text
2353 //================================================================================
2355 std::string text(netgen::NgException& ex)
2357 SMESH_Comment str("NgException");
2358 if ( strlen( netgen::multithread.task ) > 0 )
2359 str << " at " << netgen::multithread.task;
2360 str << ": " << ex.What();
2364 //================================================================================
2366 * \brief Looks for triangles lying on a SOLID
2368 //================================================================================
2370 bool hasBadElemOnSolid( const list<const SMDS_MeshElement*>& elems,
2371 SMESH_subMesh* solidSM )
2373 TopTools_IndexedMapOfShape solidSubs;
2374 TopExp::MapShapes( solidSM->GetSubShape(), solidSubs );
2375 SMESHDS_Mesh* mesh = solidSM->GetFather()->GetMeshDS();
2377 list<const SMDS_MeshElement*>::const_iterator e = elems.begin();
2378 for ( ; e != elems.end(); ++e )
2380 const SMDS_MeshElement* elem = *e;
2381 // if ( elem->GetType() != SMDSAbs_Face ) -- 23047
2383 int nbNodesOnSolid = 0, nbNodes = elem->NbNodes();
2384 SMDS_NodeIteratorPtr nIt = elem->nodeIterator();
2385 while ( nIt->more() )
2387 const SMDS_MeshNode* n = nIt->next();
2388 const TopoDS_Shape& s = mesh->IndexToShape( n->getshapeId() );
2389 nbNodesOnSolid += ( !s.IsNull() && solidSubs.Contains( s ));
2390 if ( nbNodesOnSolid > 2 ||
2391 nbNodesOnSolid == nbNodes)
2398 const double edgeMeshingTime = 0.001;
2399 const double faceMeshingTime = 0.019;
2400 const double edgeFaceMeshingTime = edgeMeshingTime + faceMeshingTime;
2401 const double faceOptimizTime = 0.06;
2402 const double voluMeshingTime = 0.15;
2403 const double volOptimizeTime = 0.77;
2406 //=============================================================================
2408 * Here we are going to use the NETGEN mesher
2410 //=============================================================================
2412 bool NETGENPlugin_Mesher::Compute()
2414 NETGENPlugin_NetgenLibWrapper ngLib;
2416 netgen::MeshingParameters& mparams = netgen::mparam;
2417 MESSAGE("Compute with:\n"
2418 " max size = " << mparams.maxh << "\n"
2419 " segments per edge = " << mparams.segmentsperedge);
2421 " growth rate = " << mparams.grading << "\n"
2422 " elements per radius = " << mparams.curvaturesafety << "\n"
2423 " second order = " << mparams.secondorder << "\n"
2424 " quad allowed = " << mparams.quad << "\n"
2425 " surface curvature = " << mparams.uselocalh << "\n"
2426 " fuse edges = " << netgen::merge_solids);
2428 SMESH_ComputeErrorPtr error = SMESH_ComputeError::New();
2429 SMESH_MesherHelper quadHelper( *_mesh );
2430 quadHelper.SetIsQuadratic( mparams.secondorder );
2432 static string debugFile = "/tmp/ngMesh.py"; /* to call toPython( _ngMesh, debugFile )
2433 while debugging netgen */
2434 // -------------------------
2435 // Prepare OCC geometry
2436 // -------------------------
2438 netgen::OCCGeometry occgeo;
2439 list< SMESH_subMesh* > meshedSM[3]; // for 0-2 dimensions
2440 NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
2441 PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
2444 _totalTime = edgeFaceMeshingTime;
2446 _totalTime += faceOptimizTime;
2448 _totalTime += voluMeshingTime + ( _optimize ? volOptimizeTime : 0 );
2449 double doneTime = 0;
2452 _curShapeIndex = -1;
2454 // -------------------------
2455 // Generate the mesh
2456 // -------------------------
2459 NETGENPlugin_ngMeshInfo initState; // it remembers size of ng mesh equal to size of Smesh
2461 SMESH_Comment comment;
2464 // vector of nodes in which node index == netgen ID
2465 vector< const SMDS_MeshNode* > nodeVec;
2473 // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
2474 mparams.uselocalh = false;
2475 mparams.grading = 0.8; // not limitited size growth
2477 if ( _simpleHyp->GetNumberOfSegments() )
2479 mparams.maxh = occgeo.boundingbox.Diam();
2482 mparams.maxh = _simpleHyp->GetLocalLength();
2485 if ( mparams.maxh == 0.0 )
2486 mparams.maxh = occgeo.boundingbox.Diam();
2487 if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
2488 mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
2490 // Local size on faces
2491 occgeo.face_maxh = mparams.maxh;
2493 // Let netgen create _ngMesh and calculate element size on not meshed shapes
2497 int startWith = netgen::MESHCONST_ANALYSE;
2498 int endWith = netgen::MESHCONST_ANALYSE;
2503 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2505 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2507 if(netgen::multithread.terminate)
2510 comment << text(err);
2512 catch (Standard_Failure& ex)
2514 comment << text(ex);
2516 err = 0; //- MESHCONST_ANALYSE isn't so important step
2519 ngLib.setMesh(( Ng_Mesh*) _ngMesh );
2521 _ngMesh->ClearFaceDescriptors(); // we make descriptors our-self
2525 // Pass 1D simple parameters to NETGEN
2526 // --------------------------------
2527 int nbSeg = _simpleHyp->GetNumberOfSegments();
2528 double segSize = _simpleHyp->GetLocalLength();
2529 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
2531 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
2533 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
2534 setLocalSize( e, segSize, *_ngMesh );
2537 else // if ( ! _simpleHyp )
2539 // Local size on vertices and edges
2540 // --------------------------------
2541 for(std::map<int,double>::const_iterator it=EdgeId2LocalSize.begin(); it!=EdgeId2LocalSize.end(); it++)
2543 int key = (*it).first;
2544 double hi = (*it).second;
2545 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2546 const TopoDS_Edge& e = TopoDS::Edge(shape);
2547 setLocalSize( e, hi, *_ngMesh );
2549 for(std::map<int,double>::const_iterator it=VertexId2LocalSize.begin(); it!=VertexId2LocalSize.end(); it++)
2551 int key = (*it).first;
2552 double hi = (*it).second;
2553 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2554 const TopoDS_Vertex& v = TopoDS::Vertex(shape);
2555 gp_Pnt p = BRep_Tool::Pnt(v);
2556 NETGENPlugin_Mesher::RestrictLocalSize( *_ngMesh, p.XYZ(), hi );
2558 for(map<int,double>::const_iterator it=FaceId2LocalSize.begin();
2559 it!=FaceId2LocalSize.end(); it++)
2561 int key = (*it).first;
2562 double val = (*it).second;
2563 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2564 int faceNgID = occgeo.fmap.FindIndex(shape);
2565 occgeo.SetFaceMaxH(faceNgID, val);
2566 for ( TopExp_Explorer edgeExp( shape, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
2567 setLocalSize( TopoDS::Edge( edgeExp.Current() ), val, *_ngMesh );
2571 // Precompute internal edges (issue 0020676) in order to
2572 // add mesh on them correctly (twice) to netgen mesh
2573 if ( !err && internals.hasInternalEdges() )
2575 // load internal shapes into OCCGeometry
2576 netgen::OCCGeometry intOccgeo;
2577 internals.getInternalEdges( intOccgeo.fmap, intOccgeo.emap, intOccgeo.vmap, meshedSM );
2578 intOccgeo.boundingbox = occgeo.boundingbox;
2579 intOccgeo.shape = occgeo.shape;
2580 intOccgeo.face_maxh.SetSize(intOccgeo.fmap.Extent());
2581 intOccgeo.face_maxh = netgen::mparam.maxh;
2582 netgen::Mesh *tmpNgMesh = NULL;
2586 // compute local H on internal shapes in the main mesh
2587 //OCCSetLocalMeshSize(intOccgeo, *_ngMesh); it deletes _ngMesh->localH
2589 // let netgen create a temporary mesh
2591 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2593 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2595 if(netgen::multithread.terminate)
2598 // copy LocalH from the main to temporary mesh
2599 initState.transferLocalH( _ngMesh, tmpNgMesh );
2601 // compute mesh on internal edges
2602 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2604 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2606 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2608 comment << text(err);
2610 catch (Standard_Failure& ex)
2612 comment << text(ex);
2615 initState.restoreLocalH( tmpNgMesh );
2617 // fill SMESH by netgen mesh
2618 vector< const SMDS_MeshNode* > tmpNodeVec;
2619 FillSMesh( intOccgeo, *tmpNgMesh, initState, *_mesh, tmpNodeVec, comment );
2620 err = ( err || !comment.empty() );
2622 nglib::Ng_DeleteMesh((nglib::Ng_Mesh*)tmpNgMesh);
2625 // Fill _ngMesh with nodes and segments of computed submeshes
2628 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_0D ]) &&
2629 FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_1D ], &quadHelper));
2631 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2636 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2641 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2643 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2645 if(netgen::multithread.terminate)
2648 comment << text(err);
2650 catch (Standard_Failure& ex)
2652 comment << text(ex);
2657 _ticTime = ( doneTime += edgeMeshingTime ) / _totalTime / _progressTic;
2659 mparams.uselocalh = true; // restore as it is used at surface optimization
2661 // ---------------------
2662 // compute surface mesh
2663 // ---------------------
2666 // Pass 2D simple parameters to NETGEN
2668 if ( double area = _simpleHyp->GetMaxElementArea() ) {
2670 mparams.maxh = sqrt(2. * area/sqrt(3.0));
2671 mparams.grading = 0.4; // moderate size growth
2674 // length from edges
2675 if ( _ngMesh->GetNSeg() ) {
2676 double edgeLength = 0;
2677 TopTools_MapOfShape visitedEdges;
2678 for ( TopExp_Explorer exp( _shape, TopAbs_EDGE ); exp.More(); exp.Next() )
2679 if( visitedEdges.Add(exp.Current()) )
2680 edgeLength += SMESH_Algo::EdgeLength( TopoDS::Edge( exp.Current() ));
2681 // we have to multiply length by 2 since for each TopoDS_Edge there
2682 // are double set of NETGEN edges, in other words, we have to
2683 // divide _ngMesh->GetNSeg() by 2.
2684 mparams.maxh = 2*edgeLength / _ngMesh->GetNSeg();
2687 mparams.maxh = 1000;
2689 mparams.grading = 0.2; // slow size growth
2691 mparams.quad = _simpleHyp->GetAllowQuadrangles();
2692 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2693 _ngMesh->SetGlobalH (mparams.maxh);
2694 netgen::Box<3> bb = occgeo.GetBoundingBox();
2695 bb.Increase (bb.Diam()/20);
2696 _ngMesh->SetLocalH (bb.PMin(), bb.PMax(), mparams.grading);
2699 // Care of vertices internal in faces (issue 0020676)
2700 if ( internals.hasInternalVertexInFace() )
2702 // store computed segments in SMESH in order not to create SMESH
2703 // edges for ng segments added by AddIntVerticesInFaces()
2704 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2705 // add segments to faces with internal vertices
2706 AddIntVerticesInFaces( occgeo, *_ngMesh, nodeVec, internals );
2707 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2710 // Build viscous layers
2711 if ( _isViscousLayers2D )
2713 if ( !internals.hasInternalVertexInFace() ) {
2714 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2715 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2717 SMESH_ProxyMesh::Ptr viscousMesh;
2718 SMESH_MesherHelper helper( *_mesh );
2719 for ( int faceID = 1; faceID <= occgeo.fmap.Extent(); ++faceID )
2721 const TopoDS_Face& F = TopoDS::Face( occgeo.fmap( faceID ));
2722 viscousMesh = StdMeshers_ViscousLayers2D::Compute( *_mesh, F );
2725 // exclude from computation ng segments built on EDGEs of F
2726 for (int i = 1; i <= _ngMesh->GetNSeg(); i++)
2728 netgen::Segment & seg = _ngMesh->LineSegment(i);
2729 if (seg.si == faceID)
2732 // add new segments to _ngMesh instead of excluded ones
2733 helper.SetSubShape( F );
2735 StdMeshers_FaceSide::GetFaceWires( F, *_mesh, /*skipMediumNodes=*/true,
2736 error, viscousMesh );
2737 error = AddSegmentsToMesh( *_ngMesh, occgeo, wires, helper, nodeVec );
2739 if ( !error ) error = SMESH_ComputeError::New();
2741 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2744 // Let netgen compute 2D mesh
2745 startWith = netgen::MESHCONST_MESHSURFACE;
2746 endWith = _optimize ? netgen::MESHCONST_OPTSURFACE : netgen::MESHCONST_MESHSURFACE;
2751 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2753 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2755 if(netgen::multithread.terminate)
2758 comment << text (err);
2760 catch (Standard_Failure& ex)
2762 comment << text(ex);
2763 //err = 1; -- try to make volumes anyway
2765 catch (netgen::NgException exc)
2767 comment << text(exc);
2768 //err = 1; -- try to make volumes anyway
2773 doneTime += faceMeshingTime + ( _optimize ? faceOptimizTime : 0 );
2774 _ticTime = doneTime / _totalTime / _progressTic;
2776 // ---------------------
2777 // generate volume mesh
2778 // ---------------------
2779 // Fill _ngMesh with nodes and faces of computed 2D submeshes
2780 if ( !err && _isVolume && ( !meshedSM[ MeshDim_2D ].empty() || mparams.quad ))
2782 // load SMESH with computed segments and faces
2783 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2785 // compute pyramids on quadrangles
2786 SMESH_ProxyMesh::Ptr proxyMesh;
2787 if ( _mesh->NbQuadrangles() > 0 )
2788 for ( int iS = 1; iS <= occgeo.somap.Extent(); ++iS )
2790 StdMeshers_QuadToTriaAdaptor* Adaptor = new StdMeshers_QuadToTriaAdaptor;
2791 proxyMesh.reset( Adaptor );
2793 int nbPyrams = _mesh->NbPyramids();
2794 Adaptor->Compute( *_mesh, occgeo.somap(iS) );
2795 if ( nbPyrams != _mesh->NbPyramids() )
2797 list< SMESH_subMesh* > quadFaceSM;
2798 for (TopExp_Explorer face(occgeo.somap(iS), TopAbs_FACE); face.More(); face.Next())
2799 if ( Adaptor->GetProxySubMesh( face.Current() ))
2801 quadFaceSM.push_back( _mesh->GetSubMesh( face.Current() ));
2802 meshedSM[ MeshDim_2D ].remove( quadFaceSM.back() );
2804 FillNgMesh(occgeo, *_ngMesh, nodeVec, quadFaceSM, &quadHelper, proxyMesh);
2807 // fill _ngMesh with faces of sub-meshes
2808 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_2D ], &quadHelper));
2809 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2810 //toPython( _ngMesh, "/tmp/ngPython.py");
2812 if (!err && _isVolume)
2814 // Pass 3D simple parameters to NETGEN
2815 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
2816 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
2818 if ( double vol = simple3d->GetMaxElementVolume() ) {
2820 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
2821 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2824 // length from faces
2825 mparams.maxh = _ngMesh->AverageH();
2827 _ngMesh->SetGlobalH (mparams.maxh);
2828 mparams.grading = 0.4;
2830 _ngMesh->CalcLocalH(mparams.grading);
2832 _ngMesh->CalcLocalH();
2835 // Care of vertices internal in solids and internal faces (issue 0020676)
2836 if ( internals.hasInternalVertexInSolid() || internals.hasInternalFaces() )
2838 // store computed faces in SMESH in order not to create SMESH
2839 // faces for ng faces added here
2840 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2841 // add ng faces to solids with internal vertices
2842 AddIntVerticesInSolids( occgeo, *_ngMesh, nodeVec, internals );
2843 // duplicate mesh faces on internal faces
2844 FixIntFaces( occgeo, *_ngMesh, internals );
2845 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2847 // Let netgen compute 3D mesh
2848 startWith = endWith = netgen::MESHCONST_MESHVOLUME;
2853 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2855 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2857 if(netgen::multithread.terminate)
2860 if ( comment.empty() ) // do not overwrite a previos error
2861 comment << text(err);
2863 catch (Standard_Failure& ex)
2865 if ( comment.empty() ) // do not overwrite a previos error
2866 comment << text(ex);
2869 catch (netgen::NgException exc)
2871 if ( comment.empty() ) // do not overwrite a previos error
2872 comment << text(exc);
2875 _ticTime = ( doneTime += voluMeshingTime ) / _totalTime / _progressTic;
2877 // Let netgen optimize 3D mesh
2878 if ( !err && _optimize )
2880 startWith = endWith = netgen::MESHCONST_OPTVOLUME;
2885 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2887 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2889 if(netgen::multithread.terminate)
2892 if ( comment.empty() ) // do not overwrite a previos error
2893 comment << text(err);
2895 catch (Standard_Failure& ex)
2897 if ( comment.empty() ) // do not overwrite a previos error
2898 comment << text(ex);
2900 catch (netgen::NgException exc)
2902 if ( comment.empty() ) // do not overwrite a previos error
2903 comment << text(exc);
2907 if (!err && mparams.secondorder > 0)
2912 if ( !meshedSM[ MeshDim_1D ].empty() )
2914 // remove segments not attached to geometry (IPAL0052479)
2915 for (int i = 1; i <= _ngMesh->GetNSeg(); ++i)
2917 const netgen::Segment & seg = _ngMesh->LineSegment (i);
2918 if ( seg.epgeominfo[ 0 ].edgenr == 0 )
2919 _ngMesh->DeleteSegment( i );
2921 _ngMesh->Compress();
2923 // convert to quadratic
2924 netgen::OCCRefinementSurfaces ref (occgeo);
2925 ref.MakeSecondOrder (*_ngMesh);
2927 // care of elements already loaded to SMESH
2928 // if ( initState._nbSegments > 0 )
2929 // makeQuadratic( occgeo.emap, _mesh );
2930 // if ( initState._nbFaces > 0 )
2931 // makeQuadratic( occgeo.fmap, _mesh );
2933 catch (Standard_Failure& ex)
2935 if ( comment.empty() ) // do not overwrite a previos error
2936 comment << "Exception in netgen at passing to 2nd order ";
2938 catch (netgen::NgException exc)
2940 if ( comment.empty() ) // do not overwrite a previos error
2941 comment << exc.What();
2946 _ticTime = 0.98 / _progressTic;
2948 int nbNod = _ngMesh->GetNP();
2949 int nbSeg = _ngMesh->GetNSeg();
2950 int nbFac = _ngMesh->GetNSE();
2951 int nbVol = _ngMesh->GetNE();
2952 bool isOK = ( !err && (_isVolume ? (nbVol > 0) : (nbFac > 0)) );
2954 MESSAGE((err ? "Mesh Generation failure" : "End of Mesh Generation") <<
2955 ", nb nodes: " << nbNod <<
2956 ", nb segments: " << nbSeg <<
2957 ", nb faces: " << nbFac <<
2958 ", nb volumes: " << nbVol);
2960 // Feed back the SMESHDS with the generated Nodes and Elements
2961 if ( true /*isOK*/ ) // get whatever built
2963 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2965 if ( quadHelper.GetIsQuadratic() ) // remove free nodes
2966 for ( size_t i = 0; i < nodeVec.size(); ++i )
2967 if ( nodeVec[i] && nodeVec[i]->NbInverseElements() == 0 )
2968 _mesh->GetMeshDS()->RemoveFreeNode( nodeVec[i], 0, /*fromGroups=*/false );
2970 SMESH_ComputeErrorPtr readErr = ReadErrors(nodeVec);
2971 if ( readErr && !readErr->myBadElements.empty() )
2974 if ( !comment.empty() && !readErr->myComment.empty() ) comment += "\n";
2975 comment += readErr->myComment;
2977 if ( error->IsOK() && ( !isOK || comment.size() > 0 ))
2978 error->myName = COMPERR_ALGO_FAILED;
2979 if ( !comment.empty() )
2980 error->myComment = comment;
2982 // SetIsAlwaysComputed( true ) to empty sub-meshes, which
2983 // appear if the geometry contains coincident sub-shape due
2984 // to bool merge_solids = 1; in netgen/libsrc/occ/occgenmesh.cpp
2985 const int nbMaps = 2;
2986 const TopTools_IndexedMapOfShape* geoMaps[nbMaps] =
2987 { & occgeo.vmap, & occgeo.emap/*, & occgeo.fmap*/ };
2988 for ( int iMap = 0; iMap < nbMaps; ++iMap )
2989 for (int i = 1; i <= geoMaps[iMap]->Extent(); i++)
2990 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( geoMaps[iMap]->FindKey(i)))
2991 if ( !sm->IsMeshComputed() )
2992 sm->SetIsAlwaysComputed( true );
2994 // set bad compute error to subshapes of all failed sub-shapes
2995 if ( !error->IsOK() )
2997 bool pb2D = false, pb3D = false;
2998 for (int i = 1; i <= occgeo.fmap.Extent(); i++) {
2999 int status = occgeo.facemeshstatus[i-1];
3000 if (status == netgen::FACE_MESHED_OK ) continue;
3001 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.fmap( i ))) {
3002 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
3003 if ( !smError || smError->IsOK() ) {
3004 if ( status == netgen::FACE_FAILED )
3005 smError.reset( new SMESH_ComputeError( *error ));
3007 smError.reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED, "Ignored" ));
3008 if ( SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
3009 smError->myName = COMPERR_WARNING;
3011 pb2D = pb2D || smError->IsKO();
3014 if ( !pb2D ) // all faces are OK
3015 for (int i = 1; i <= occgeo.somap.Extent(); i++)
3016 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.somap( i )))
3018 bool smComputed = nbVol && !sm->IsEmpty();
3019 if ( smComputed && internals.hasInternalVertexInSolid( sm->GetId() ))
3021 int nbIntV = internals.getSolidsWithVertices().find( sm->GetId() )->second.size();
3022 SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
3023 smComputed = ( smDS->NbElements() > 0 || smDS->NbNodes() > nbIntV );
3025 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
3026 if ( !smComputed && ( !smError || smError->IsOK() ))
3028 smError.reset( new SMESH_ComputeError( *error ));
3029 if ( nbVol && SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
3031 smError->myName = COMPERR_WARNING;
3033 else if ( !smError->myBadElements.empty() ) // bad surface mesh
3035 if ( !hasBadElemOnSolid( smError->myBadElements, sm ))
3039 pb3D = pb3D || ( smError && smError->IsKO() );
3041 if ( !pb2D && !pb3D )
3042 err = 0; // no fatal errors, only warnings
3045 ngLib._isComputeOk = !err;
3050 //=============================================================================
3054 //=============================================================================
3055 bool NETGENPlugin_Mesher::Evaluate(MapShapeNbElems& aResMap)
3057 netgen::MeshingParameters& mparams = netgen::mparam;
3060 // -------------------------
3061 // Prepare OCC geometry
3062 // -------------------------
3063 netgen::OCCGeometry occgeo;
3064 list< SMESH_subMesh* > meshedSM[4]; // for 0-3 dimensions
3065 NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
3066 PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
3068 bool tooManyElems = false;
3069 const int hugeNb = std::numeric_limits<int>::max() / 100;
3074 // pass 1D simple parameters to NETGEN
3077 // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
3078 mparams.uselocalh = false;
3079 mparams.grading = 0.8; // not limitited size growth
3081 if ( _simpleHyp->GetNumberOfSegments() )
3083 mparams.maxh = occgeo.boundingbox.Diam();
3086 mparams.maxh = _simpleHyp->GetLocalLength();
3089 if ( mparams.maxh == 0.0 )
3090 mparams.maxh = occgeo.boundingbox.Diam();
3091 if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
3092 mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
3094 // let netgen create _ngMesh and calculate element size on not meshed shapes
3095 NETGENPlugin_NetgenLibWrapper ngLib;
3096 netgen::Mesh *ngMesh = NULL;
3100 int startWith = netgen::MESHCONST_ANALYSE;
3101 int endWith = netgen::MESHCONST_MESHEDGES;
3103 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, mparams, startWith, endWith);
3105 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, startWith, endWith, optstr);
3108 if(netgen::multithread.terminate)
3111 ngLib.setMesh(( Ng_Mesh*) ngMesh );
3113 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( _shape ))
3114 sm->GetComputeError().reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED ));
3119 // Pass 1D simple parameters to NETGEN
3120 // --------------------------------
3121 int nbSeg = _simpleHyp->GetNumberOfSegments();
3122 double segSize = _simpleHyp->GetLocalLength();
3123 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
3125 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
3127 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
3128 setLocalSize( e, segSize, *ngMesh );
3131 else // if ( ! _simpleHyp )
3133 // Local size on vertices and edges
3134 // --------------------------------
3135 for(std::map<int,double>::const_iterator it=EdgeId2LocalSize.begin(); it!=EdgeId2LocalSize.end(); it++)
3137 int key = (*it).first;
3138 double hi = (*it).second;
3139 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3140 const TopoDS_Edge& e = TopoDS::Edge(shape);
3141 setLocalSize( e, hi, *ngMesh );
3143 for(std::map<int,double>::const_iterator it=VertexId2LocalSize.begin(); it!=VertexId2LocalSize.end(); it++)
3145 int key = (*it).first;
3146 double hi = (*it).second;
3147 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3148 const TopoDS_Vertex& v = TopoDS::Vertex(shape);
3149 gp_Pnt p = BRep_Tool::Pnt(v);
3150 NETGENPlugin_Mesher::RestrictLocalSize( *ngMesh, p.XYZ(), hi );
3152 for(map<int,double>::const_iterator it=FaceId2LocalSize.begin();
3153 it!=FaceId2LocalSize.end(); it++)
3155 int key = (*it).first;
3156 double val = (*it).second;
3157 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3158 int faceNgID = occgeo.fmap.FindIndex(shape);
3159 occgeo.SetFaceMaxH(faceNgID, val);
3160 for ( TopExp_Explorer edgeExp( shape, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
3161 setLocalSize( TopoDS::Edge( edgeExp.Current() ), val, *ngMesh );
3164 // calculate total nb of segments and length of edges
3165 double fullLen = 0.0;
3167 int entity = mparams.secondorder > 0 ? SMDSEntity_Quad_Edge : SMDSEntity_Edge;
3168 TopTools_DataMapOfShapeInteger Edge2NbSeg;
3169 for (TopExp_Explorer exp(_shape, TopAbs_EDGE); exp.More(); exp.Next())
3171 TopoDS_Edge E = TopoDS::Edge( exp.Current() );
3172 if( !Edge2NbSeg.Bind(E,0) )
3175 double aLen = SMESH_Algo::EdgeLength(E);
3178 vector<int>& aVec = aResMap[_mesh->GetSubMesh(E)];
3180 aVec.resize( SMDSEntity_Last, 0);
3182 fullNbSeg += aVec[ entity ];
3185 // store nb of segments computed by Netgen
3186 NCollection_Map<Link> linkMap;
3187 for (int i = 1; i <= ngMesh->GetNSeg(); ++i )
3189 const netgen::Segment& seg = ngMesh->LineSegment(i);
3190 Link link(seg[0], seg[1]);
3191 if ( !linkMap.Add( link )) continue;
3192 int aGeomEdgeInd = seg.epgeominfo[0].edgenr;
3193 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
3195 vector<int>& aVec = aResMap[_mesh->GetSubMesh(occgeo.emap(aGeomEdgeInd))];
3199 // store nb of nodes on edges computed by Netgen
3200 TopTools_DataMapIteratorOfDataMapOfShapeInteger Edge2NbSegIt(Edge2NbSeg);
3201 for (; Edge2NbSegIt.More(); Edge2NbSegIt.Next())
3203 vector<int>& aVec = aResMap[_mesh->GetSubMesh(Edge2NbSegIt.Key())];
3204 if ( aVec[ entity ] > 1 && aVec[ SMDSEntity_Node ] == 0 )
3205 aVec[SMDSEntity_Node] = mparams.secondorder > 0 ? 2*aVec[ entity ]-1 : aVec[ entity ]-1;
3207 fullNbSeg += aVec[ entity ];
3208 Edge2NbSeg( Edge2NbSegIt.Key() ) = aVec[ entity ];
3210 if ( fullNbSeg == 0 )
3217 if ( double area = _simpleHyp->GetMaxElementArea() ) {
3219 mparams.maxh = sqrt(2. * area/sqrt(3.0));
3220 mparams.grading = 0.4; // moderate size growth
3223 // length from edges
3224 mparams.maxh = fullLen/fullNbSeg;
3225 mparams.grading = 0.2; // slow size growth
3228 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3229 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3231 for (TopExp_Explorer exp(_shape, TopAbs_FACE); exp.More(); exp.Next())
3233 TopoDS_Face F = TopoDS::Face( exp.Current() );
3234 SMESH_subMesh *sm = _mesh->GetSubMesh(F);
3236 BRepGProp::SurfaceProperties(F,G);
3237 double anArea = G.Mass();
3238 tooManyElems = tooManyElems || ( anArea/hugeNb > mparams.maxh*mparams.maxh );
3240 if ( !tooManyElems )
3242 TopTools_MapOfShape egdes;
3243 for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next())
3244 if ( egdes.Add( exp1.Current() ))
3245 nb1d += Edge2NbSeg.Find(exp1.Current());
3247 int nbFaces = tooManyElems ? hugeNb : int( 4*anArea / (mparams.maxh*mparams.maxh*sqrt(3.)));
3248 int nbNodes = tooManyElems ? hugeNb : (( nbFaces*3 - (nb1d-1)*2 ) / 6 + 1 );
3250 vector<int> aVec(SMDSEntity_Last, 0);
3251 if( mparams.secondorder > 0 ) {
3252 int nb1d_in = (nbFaces*3 - nb1d) / 2;
3253 aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
3254 aVec[SMDSEntity_Quad_Triangle] = nbFaces;
3257 aVec[SMDSEntity_Node] = Max ( nbNodes, 0 );
3258 aVec[SMDSEntity_Triangle] = nbFaces;
3260 aResMap[sm].swap(aVec);
3267 // pass 3D simple parameters to NETGEN
3268 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
3269 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
3271 if ( double vol = simple3d->GetMaxElementVolume() ) {
3273 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
3274 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3277 // using previous length from faces
3279 mparams.grading = 0.4;
3280 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3283 BRepGProp::VolumeProperties(_shape,G);
3284 double aVolume = G.Mass();
3285 double tetrVol = 0.1179*mparams.maxh*mparams.maxh*mparams.maxh;
3286 tooManyElems = tooManyElems || ( aVolume/hugeNb > tetrVol );
3287 int nbVols = tooManyElems ? hugeNb : int(aVolume/tetrVol);
3288 int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
3289 vector<int> aVec(SMDSEntity_Last, 0 );
3290 if ( tooManyElems ) // avoid FPE
3292 aVec[SMDSEntity_Node] = hugeNb;
3293 aVec[ mparams.secondorder > 0 ? SMDSEntity_Quad_Tetra : SMDSEntity_Tetra] = hugeNb;
3297 if( mparams.secondorder > 0 ) {
3298 aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
3299 aVec[SMDSEntity_Quad_Tetra] = nbVols;
3302 aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
3303 aVec[SMDSEntity_Tetra] = nbVols;
3306 SMESH_subMesh *sm = _mesh->GetSubMesh(_shape);
3307 aResMap[sm].swap(aVec);
3313 double NETGENPlugin_Mesher::GetProgress(const SMESH_Algo* holder,
3314 const int * algoProgressTic,
3315 const double * algoProgress) const
3317 ((int&) _progressTic ) = *algoProgressTic + 1;
3319 if ( !_occgeom ) return 0;
3321 double progress = -1;
3324 if ( _ticTime < 0 && netgen::multithread.task[0] == 'O'/*Optimizing surface*/ )
3326 ((double&) _ticTime ) = edgeFaceMeshingTime / _totalTime / _progressTic;
3328 else if ( !_optimize /*&& _occgeom->fmap.Extent() > 1*/ )
3330 int doneShapeIndex = -1;
3331 while ( doneShapeIndex+1 < _occgeom->facemeshstatus.Size() &&
3332 _occgeom->facemeshstatus[ doneShapeIndex+1 ])
3334 if ( doneShapeIndex+1 != _curShapeIndex )
3336 ((int&) _curShapeIndex) = doneShapeIndex+1;
3337 double doneShapeRate = _curShapeIndex / double( _occgeom->fmap.Extent() );
3338 double doneTime = edgeMeshingTime + doneShapeRate * faceMeshingTime;
3339 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3340 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3341 // << " " << doneTime / _totalTime / _progressTic << endl;
3345 else if ( !_optimize && _occgeom->somap.Extent() > 1 )
3347 int curShapeIndex = _curShapeIndex;
3348 if ( _ngMesh->GetNE() > 0 )
3350 netgen::Element el = (*_ngMesh)[netgen::ElementIndex( _ngMesh->GetNE()-1 )];
3351 curShapeIndex = el.GetIndex();
3353 if ( curShapeIndex != _curShapeIndex )
3355 ((int&) _curShapeIndex) = curShapeIndex;
3356 double doneShapeRate = _curShapeIndex / double( _occgeom->somap.Extent() );
3357 double doneTime = edgeFaceMeshingTime + doneShapeRate * voluMeshingTime;
3358 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3359 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3360 // << " " << doneTime / _totalTime / _progressTic << endl;
3364 progress = Max( *algoProgressTic * _ticTime, *algoProgress );
3367 ((int&) *algoProgressTic )++;
3368 ((double&) *algoProgress) = progress;
3370 //cout << progress << " " << *algoProgressTic << " " << netgen::multithread.task << " "<< _ticTime << endl;
3372 return Min( progress, 0.99 );
3375 //================================================================================
3377 * \brief Remove "test.out" and "problemfaces" files in current directory
3379 //================================================================================
3381 void NETGENPlugin_Mesher::RemoveTmpFiles()
3383 bool rm = SMESH_File("test.out").remove() ;
3385 if (rm && netgen::testout)
3387 delete netgen::testout;
3388 netgen::testout = 0;
3391 SMESH_File("problemfaces").remove();
3392 SMESH_File("occmesh.rep").remove();
3395 //================================================================================
3397 * \brief Read mesh entities preventing successful computation from "test.out" file
3399 //================================================================================
3401 SMESH_ComputeErrorPtr
3402 NETGENPlugin_Mesher::ReadErrors(const vector<const SMDS_MeshNode* >& nodeVec)
3404 SMESH_ComputeErrorPtr err = SMESH_ComputeError::New
3405 (COMPERR_BAD_INPUT_MESH, "Some edges multiple times in surface mesh");
3406 SMESH_File file("test.out");
3408 vector<int> three1(3), three2(3);
3409 const char* badEdgeStr = " multiple times in surface mesh";
3410 const int badEdgeStrLen = strlen( badEdgeStr );
3411 const int nbNodes = nodeVec.size();
3413 while( !file.eof() )
3415 if ( strncmp( file, "Edge ", 5 ) == 0 &&
3416 file.getInts( two ) &&
3417 strncmp( file, badEdgeStr, badEdgeStrLen ) == 0 &&
3418 two[0] < nbNodes && two[1] < nbNodes )
3420 err->myBadElements.push_back( new SMDS_LinearEdge( nodeVec[ two[0]], nodeVec[ two[1]] ));
3421 file += badEdgeStrLen;
3423 else if ( strncmp( file, "Intersecting: ", 14 ) == 0 )
3426 // openelement 18 with open element 126
3430 const char* pos = file;
3431 bool ok = ( strncmp( file, "openelement ", 12 ) == 0 );
3432 ok = ok && file.getInts( two );
3433 ok = ok && file.getInts( three1 );
3434 ok = ok && file.getInts( three2 );
3435 for ( int i = 0; ok && i < 3; ++i )
3436 ok = ( three1[i] < nbNodes && nodeVec[ three1[i]]);
3437 for ( int i = 0; ok && i < 3; ++i )
3438 ok = ( three2[i] < nbNodes && nodeVec[ three2[i]]);
3441 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three1[0]],
3442 nodeVec[ three1[1]],
3443 nodeVec[ three1[2]]));
3444 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three2[0]],
3445 nodeVec[ three2[1]],
3446 nodeVec[ three2[2]]));
3447 err->myComment = "Intersecting triangles";
3461 size_t nbBadElems = err->myBadElements.size();
3468 //================================================================================
3470 * \brief Write a python script creating an equivalent SALOME mesh.
3471 * This is useful to see what mesh is passed as input for the next step of mesh
3472 * generation (of mesh of higher dimension)
3474 //================================================================================
3476 void NETGENPlugin_Mesher::toPython( const netgen::Mesh* ngMesh,
3477 const std::string& pyFile)
3479 ofstream outfile(pyFile.c_str(), ios::out);
3480 if ( !outfile ) return;
3482 outfile << "import SMESH" << endl
3483 << "from salome.smesh import smeshBuilder" << endl
3484 << "smesh = smeshBuilder.New(salome.myStudy)" << endl
3485 << "mesh = smesh.Mesh()" << endl << endl;
3487 using namespace netgen;
3489 for (pi = PointIndex::BASE;
3490 pi < ngMesh->GetNP()+PointIndex::BASE; pi++)
3492 outfile << "mesh.AddNode( ";
3493 outfile << (*ngMesh)[pi](0) << ", ";
3494 outfile << (*ngMesh)[pi](1) << ", ";
3495 outfile << (*ngMesh)[pi](2) << ") ## "<< pi << endl;
3498 int nbDom = ngMesh->GetNDomains();
3499 for ( int i = 0; i < nbDom; ++i )
3500 outfile<< "grp" << i+1 << " = mesh.CreateEmptyGroup( SMESH.FACE, 'domain"<< i+1 << "')"<< endl;
3502 SurfaceElementIndex sei;
3503 for (sei = 0; sei < ngMesh->GetNSE(); sei++)
3505 outfile << "mesh.AddFace([ ";
3506 Element2d sel = (*ngMesh)[sei];
3507 for (int j = 0; j < sel.GetNP(); j++)
3508 outfile << sel[j] << ( j+1 < sel.GetNP() ? ", " : " ])");
3509 if ( sel.IsDeleted() ) outfile << " ## IsDeleted ";
3512 if ((*ngMesh)[sei].GetIndex())
3514 if ( int dom1 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainIn())
3515 outfile << "grp"<< dom1 <<".Add([ " << (int)sei+1 << " ])" << endl;
3516 if ( int dom2 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainOut())
3517 outfile << "grp"<< dom2 <<".Add([ " << (int)sei+1 << " ])" << endl;
3521 for (ElementIndex ei = 0; ei < ngMesh->GetNE(); ei++)
3523 Element el = (*ngMesh)[ei];
3524 outfile << "mesh.AddVolume([ ";
3525 for (int j = 0; j < el.GetNP(); j++)
3526 outfile << el[j] << ( j+1 < el.GetNP() ? ", " : " ])");
3530 for (int i = 1; i <= ngMesh->GetNSeg(); i++)
3532 const Segment & seg = ngMesh->LineSegment (i);
3533 outfile << "mesh.AddEdge([ "
3535 << seg[1] << " ])" << endl;
3537 cout << "Write " << pyFile << endl;
3540 //================================================================================
3542 * \brief Constructor of NETGENPlugin_ngMeshInfo
3544 //================================================================================
3546 NETGENPlugin_ngMeshInfo::NETGENPlugin_ngMeshInfo( netgen::Mesh* ngMesh):
3551 _nbNodes = ngMesh->GetNP();
3552 _nbSegments = ngMesh->GetNSeg();
3553 _nbFaces = ngMesh->GetNSE();
3554 _nbVolumes = ngMesh->GetNE();
3558 _nbNodes = _nbSegments = _nbFaces = _nbVolumes = 0;
3562 //================================================================================
3564 * \brief Copy LocalH member from one netgen mesh to another
3566 //================================================================================
3568 void NETGENPlugin_ngMeshInfo::transferLocalH( netgen::Mesh* fromMesh,
3569 netgen::Mesh* toMesh )
3571 if ( !fromMesh->LocalHFunctionGenerated() ) return;
3572 if ( !toMesh->LocalHFunctionGenerated() )
3574 toMesh->CalcLocalH(netgen::mparam.grading);
3576 toMesh->CalcLocalH();
3579 const size_t size = sizeof( netgen::LocalH );
3580 _copyOfLocalH = new char[ size ];
3581 memcpy( (void*)_copyOfLocalH, (void*)&toMesh->LocalHFunction(), size );
3582 memcpy( (void*)&toMesh->LocalHFunction(), (void*)&fromMesh->LocalHFunction(), size );
3585 //================================================================================
3587 * \brief Restore LocalH member of a netgen mesh
3589 //================================================================================
3591 void NETGENPlugin_ngMeshInfo::restoreLocalH( netgen::Mesh* toMesh )
3593 if ( _copyOfLocalH )
3595 const size_t size = sizeof( netgen::LocalH );
3596 memcpy( (void*)&toMesh->LocalHFunction(), (void*)_copyOfLocalH, size );
3597 delete [] _copyOfLocalH;
3602 //================================================================================
3604 * \brief Find "internal" sub-shapes
3606 //================================================================================
3608 NETGENPlugin_Internals::NETGENPlugin_Internals( SMESH_Mesh& mesh,
3609 const TopoDS_Shape& shape,
3611 : _mesh( mesh ), _is3D( is3D )
3613 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3615 TopExp_Explorer f,e;
3616 for ( f.Init( shape, TopAbs_FACE ); f.More(); f.Next() )
3618 int faceID = meshDS->ShapeToIndex( f.Current() );
3620 // find not computed internal edges
3622 for ( e.Init( f.Current().Oriented(TopAbs_FORWARD), TopAbs_EDGE ); e.More(); e.Next() )
3623 if ( e.Current().Orientation() == TopAbs_INTERNAL )
3625 SMESH_subMesh* eSM = mesh.GetSubMesh( e.Current() );
3626 if ( eSM->IsEmpty() )
3628 _e2face.insert( make_pair( eSM->GetId(), faceID ));
3629 for ( TopoDS_Iterator v(e.Current()); v.More(); v.Next() )
3630 _e2face.insert( make_pair( meshDS->ShapeToIndex( v.Value() ), faceID ));
3634 // find internal vertices in a face
3635 set<int> intVV; // issue 0020850 where same vertex is twice in a face
3636 for ( TopoDS_Iterator fSub( f.Current() ); fSub.More(); fSub.Next())
3637 if ( fSub.Value().ShapeType() == TopAbs_VERTEX )
3639 int vID = meshDS->ShapeToIndex( fSub.Value() );
3640 if ( intVV.insert( vID ).second )
3641 _f2v[ faceID ].push_back( vID );
3646 // find internal faces and their subshapes where nodes are to be doubled
3647 // to make a crack with non-sewed borders
3649 if ( f.Current().Orientation() == TopAbs_INTERNAL )
3651 _intShapes.insert( meshDS->ShapeToIndex( f.Current() ));
3654 list< TopoDS_Shape > edges;
3655 for ( e.Init( f.Current(), TopAbs_EDGE ); e.More(); e.Next())
3656 if ( SMESH_MesherHelper::NbAncestors( e.Current(), mesh, TopAbs_FACE ) > 1 )
3658 _intShapes.insert( meshDS->ShapeToIndex( e.Current() ));
3659 edges.push_back( e.Current() );
3660 // find border faces
3661 PShapeIteratorPtr fIt =
3662 SMESH_MesherHelper::GetAncestors( edges.back(),mesh,TopAbs_FACE );
3663 while ( const TopoDS_Shape* pFace = fIt->next() )
3664 if ( !pFace->IsSame( f.Current() ))
3665 _borderFaces.insert( meshDS->ShapeToIndex( *pFace ));
3668 // we consider vertex internal if it is shared by more than one internal edge
3669 list< TopoDS_Shape >::iterator edge = edges.begin();
3670 for ( ; edge != edges.end(); ++edge )
3671 for ( TopoDS_Iterator v( *edge ); v.More(); v.Next() )
3673 set<int> internalEdges;
3674 PShapeIteratorPtr eIt =
3675 SMESH_MesherHelper::GetAncestors( v.Value(),mesh,TopAbs_EDGE );
3676 while ( const TopoDS_Shape* pEdge = eIt->next() )
3678 int edgeID = meshDS->ShapeToIndex( *pEdge );
3679 if ( isInternalShape( edgeID ))
3680 internalEdges.insert( edgeID );
3682 if ( internalEdges.size() > 1 )
3683 _intShapes.insert( meshDS->ShapeToIndex( v.Value() ));
3687 } // loop on geom faces
3689 // find vertices internal in solids
3692 for ( TopExp_Explorer so(shape, TopAbs_SOLID); so.More(); so.Next())
3694 int soID = meshDS->ShapeToIndex( so.Current() );
3695 for ( TopoDS_Iterator soSub( so.Current() ); soSub.More(); soSub.Next())
3696 if ( soSub.Value().ShapeType() == TopAbs_VERTEX )
3697 _s2v[ soID ].push_back( meshDS->ShapeToIndex( soSub.Value() ));
3702 //================================================================================
3704 * \brief Find mesh faces on non-internal geom faces sharing internal edge
3705 * some nodes of which are to be doubled to make the second border of the "crack"
3707 //================================================================================
3709 void NETGENPlugin_Internals::findBorderElements( TIDSortedElemSet & borderElems )
3711 if ( _intShapes.empty() ) return;
3713 SMESH_Mesh& mesh = const_cast<SMESH_Mesh&>(_mesh);
3714 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3716 // loop on internal geom edges
3717 set<int>::const_iterator intShapeId = _intShapes.begin();
3718 for ( ; intShapeId != _intShapes.end(); ++intShapeId )
3720 const TopoDS_Shape& s = meshDS->IndexToShape( *intShapeId );
3721 if ( s.ShapeType() != TopAbs_EDGE ) continue;
3723 // get internal and non-internal geom faces sharing the internal edge <s>
3725 set<int>::iterator bordFace = _borderFaces.end();
3726 PShapeIteratorPtr faces = SMESH_MesherHelper::GetAncestors( s, _mesh, TopAbs_FACE );
3727 while ( const TopoDS_Shape* pFace = faces->next() )
3729 int faceID = meshDS->ShapeToIndex( *pFace );
3730 if ( isInternalShape( faceID ))
3733 bordFace = _borderFaces.insert( faceID ).first;
3735 if ( bordFace == _borderFaces.end() || !intFace ) continue;
3737 // get all links of mesh faces on internal geom face sharing nodes on edge <s>
3738 set< SMESH_OrientedLink > links; //!< links of faces on internal geom face
3739 list<const SMDS_MeshElement*> suspectFaces[2]; //!< mesh faces on border geom faces
3740 int nbSuspectFaces = 0;
3741 SMESHDS_SubMesh* intFaceSM = meshDS->MeshElements( intFace );
3742 if ( !intFaceSM || intFaceSM->NbElements() == 0 ) continue;
3743 SMESH_subMeshIteratorPtr smIt = mesh.GetSubMesh( s )->getDependsOnIterator(true,true);
3744 while ( smIt->more() )
3746 SMESHDS_SubMesh* sm = smIt->next()->GetSubMeshDS();
3747 if ( !sm ) continue;
3748 SMDS_NodeIteratorPtr nIt = sm->GetNodes();
3749 while ( nIt->more() )
3751 const SMDS_MeshNode* nOnEdge = nIt->next();
3752 SMDS_ElemIteratorPtr fIt = nOnEdge->GetInverseElementIterator(SMDSAbs_Face);
3753 while ( fIt->more() )
3755 const SMDS_MeshElement* f = fIt->next();
3756 const int nbNodes = f->NbCornerNodes();
3757 if ( intFaceSM->Contains( f ))
3759 for ( int i = 0; i < nbNodes; ++i )
3760 links.insert( SMESH_OrientedLink( f->GetNode(i), f->GetNode((i+1)%nbNodes)));
3765 for ( int i = 0; i < nbNodes; ++i )
3766 nbDblNodes += isInternalShape( f->GetNode(i)->getshapeId() );
3768 suspectFaces[ nbDblNodes < 2 ].push_back( f );
3774 // suspectFaces[0] having link with same orientation as mesh faces on
3775 // the internal geom face are <borderElems>. suspectFaces[1] have
3776 // only one node on edge <s>, we decide on them later (at the 2nd loop)
3777 // by links of <borderElems> found at the 1st and 2nd loops
3778 set< SMESH_OrientedLink > borderLinks;
3779 for ( int isPostponed = 0; isPostponed < 2; ++isPostponed )
3781 list<const SMDS_MeshElement*>::iterator fIt = suspectFaces[isPostponed].begin();
3782 for ( int nbF = 0; fIt != suspectFaces[isPostponed].end(); ++fIt, ++nbF )
3784 const SMDS_MeshElement* f = *fIt;
3785 bool isBorder = false, linkFound = false, borderLinkFound = false;
3786 list< SMESH_OrientedLink > faceLinks;
3787 int nbNodes = f->NbCornerNodes();
3788 for ( int i = 0; i < nbNodes; ++i )
3790 SMESH_OrientedLink link( f->GetNode(i), f->GetNode((i+1)%nbNodes));
3791 faceLinks.push_back( link );
3794 set< SMESH_OrientedLink >::iterator foundLink = links.find( link );
3795 if ( foundLink != links.end() )
3798 isBorder = ( foundLink->_reversed == link._reversed );
3799 if ( !isBorder && !isPostponed ) break;
3800 faceLinks.pop_back();
3802 else if ( isPostponed && !borderLinkFound )
3804 foundLink = borderLinks.find( link );
3805 if ( foundLink != borderLinks.end() )
3807 borderLinkFound = true;
3808 isBorder = ( foundLink->_reversed != link._reversed );
3815 borderElems.insert( f );
3816 borderLinks.insert( faceLinks.begin(), faceLinks.end() );
3818 else if ( !linkFound && !borderLinkFound )
3820 suspectFaces[1].push_back( f );
3821 if ( nbF > 2 * nbSuspectFaces )
3822 break; // dead loop protection
3829 //================================================================================
3831 * \brief put internal shapes in maps and fill in submeshes to precompute
3833 //================================================================================
3835 void NETGENPlugin_Internals::getInternalEdges( TopTools_IndexedMapOfShape& fmap,
3836 TopTools_IndexedMapOfShape& emap,
3837 TopTools_IndexedMapOfShape& vmap,
3838 list< SMESH_subMesh* > smToPrecompute[])
3840 if ( !hasInternalEdges() ) return;
3841 map<int,int>::const_iterator ev_face = _e2face.begin();
3842 for ( ; ev_face != _e2face.end(); ++ev_face )
3844 const TopoDS_Shape& ev = _mesh.GetMeshDS()->IndexToShape( ev_face->first );
3845 const TopoDS_Shape& face = _mesh.GetMeshDS()->IndexToShape( ev_face->second );
3847 ( ev.ShapeType() == TopAbs_EDGE ? emap : vmap ).Add( ev );
3849 //cout<<"INTERNAL EDGE or VERTEX "<<ev_face->first<<" on face "<<ev_face->second<<endl;
3851 smToPrecompute[ MeshDim_1D ].push_back( _mesh.GetSubMeshContaining( ev_face->first ));
3855 //================================================================================
3857 * \brief return shapes and submeshes to be meshed and already meshed boundary submeshes
3859 //================================================================================
3861 void NETGENPlugin_Internals::getInternalFaces( TopTools_IndexedMapOfShape& fmap,
3862 TopTools_IndexedMapOfShape& emap,
3863 list< SMESH_subMesh* >& intFaceSM,
3864 list< SMESH_subMesh* >& boundarySM)
3866 if ( !hasInternalFaces() ) return;
3868 // <fmap> and <emap> are for not yet meshed shapes
3869 // <intFaceSM> is for submeshes of faces
3870 // <boundarySM> is for meshed edges and vertices
3875 set<int> shapeIDs ( _intShapes );
3876 if ( !_borderFaces.empty() )
3877 shapeIDs.insert( _borderFaces.begin(), _borderFaces.end() );
3879 set<int>::const_iterator intS = shapeIDs.begin();
3880 for ( ; intS != shapeIDs.end(); ++intS )
3882 SMESH_subMesh* sm = _mesh.GetSubMeshContaining( *intS );
3884 if ( sm->GetSubShape().ShapeType() != TopAbs_FACE ) continue;
3886 intFaceSM.push_back( sm );
3888 // add submeshes of not computed internal faces
3889 if ( !sm->IsEmpty() ) continue;
3891 SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(true,true);
3892 while ( smIt->more() )
3895 const TopoDS_Shape& s = sm->GetSubShape();
3897 if ( sm->IsEmpty() )
3900 switch ( s.ShapeType() ) {
3901 case TopAbs_FACE: fmap.Add ( s ); break;
3902 case TopAbs_EDGE: emap.Add ( s ); break;
3908 if ( s.ShapeType() != TopAbs_FACE )
3909 boundarySM.push_back( sm );
3915 //================================================================================
3917 * \brief Return true if given shape is to be precomputed in order to be correctly
3918 * added to netgen mesh
3920 //================================================================================
3922 bool NETGENPlugin_Internals::isShapeToPrecompute(const TopoDS_Shape& s)
3924 int shapeID = _mesh.GetMeshDS()->ShapeToIndex( s );
3925 switch ( s.ShapeType() ) {
3926 case TopAbs_FACE : break; //return isInternalShape( shapeID ) || isBorderFace( shapeID );
3927 case TopAbs_EDGE : return isInternalEdge( shapeID );
3928 case TopAbs_VERTEX: break;
3934 //================================================================================
3936 * \brief Return SMESH
3938 //================================================================================
3940 SMESH_Mesh& NETGENPlugin_Internals::getMesh() const
3942 return const_cast<SMESH_Mesh&>( _mesh );
3945 //================================================================================
3947 * \brief Initialize netgen library
3949 //================================================================================
3951 NETGENPlugin_NetgenLibWrapper::NETGENPlugin_NetgenLibWrapper()
3955 _isComputeOk = false;
3957 if ( !getenv( "KEEP_NETGEN_OUTPUT" ))
3959 // redirect all netgen output (mycout,myerr,cout) to _outputFileName
3960 _outputFileName = getOutputFileName();
3961 netgen::mycout = new ofstream ( _outputFileName.c_str() );
3962 netgen::myerr = netgen::mycout;
3963 _coutBuffer = std::cout.rdbuf();
3965 cout << "NOTE: netgen output is redirected to file " << _outputFileName << endl;
3967 std::cout.rdbuf( netgen::mycout->rdbuf() );
3971 _ngMesh = Ng_NewMesh();
3974 //================================================================================
3976 * \brief Finish using netgen library
3978 //================================================================================
3980 NETGENPlugin_NetgenLibWrapper::~NETGENPlugin_NetgenLibWrapper()
3982 Ng_DeleteMesh( _ngMesh );
3984 NETGENPlugin_Mesher::RemoveTmpFiles();
3986 std::cout.rdbuf( _coutBuffer );
3993 //================================================================================
3995 * \brief Set netgen mesh to delete at destruction
3997 //================================================================================
3999 void NETGENPlugin_NetgenLibWrapper::setMesh( Ng_Mesh* mesh )
4002 Ng_DeleteMesh( _ngMesh );
4006 //================================================================================
4008 * \brief Return a unique file name
4010 //================================================================================
4012 std::string NETGENPlugin_NetgenLibWrapper::getOutputFileName()
4014 std::string aTmpDir = SALOMEDS_Tool::GetTmpDir();
4016 TCollection_AsciiString aGenericName = (char*)aTmpDir.c_str();
4017 aGenericName += "NETGEN_";
4019 aGenericName += getpid();
4021 aGenericName += _getpid();
4023 aGenericName += "_";
4024 aGenericName += Abs((Standard_Integer)(long) aGenericName.ToCString());
4025 aGenericName += ".out";
4027 return aGenericName.ToCString();
4030 //================================================================================
4032 * \brief Remove file with netgen output
4034 //================================================================================
4036 void NETGENPlugin_NetgenLibWrapper::removeOutputFile()
4038 if ( !_outputFileName.empty() )
4040 if ( netgen::mycout )
4042 delete netgen::mycout;
4046 string tmpDir = SALOMEDS_Tool::GetDirFromPath ( _outputFileName );
4047 string aFileName = SALOMEDS_Tool::GetNameFromPath( _outputFileName ) + ".out";
4048 SALOMEDS::ListOfFileNames_var aFiles = new SALOMEDS::ListOfFileNames;
4050 aFiles[0] = aFileName.c_str();
4052 SALOMEDS_Tool::RemoveTemporaryFiles( tmpDir.c_str(), aFiles.in(), true );