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 netgen::Point<3> p1 = netgen::Point<3> (x1,y1,z1);
576 netgen::Point<3> p2 = netgen::Point<3> (x2,y2,z2);
577 occgeo.boundingbox = netgen::Box<3> (p1,p2);
579 occgeo.shape = shape;
582 // fill maps of shapes of occgeo with not yet meshed subshapes
584 // get root submeshes
585 list< SMESH_subMesh* > rootSM;
586 const int shapeID = mesh.GetMeshDS()->ShapeToIndex( shape );
587 if ( shapeID > 0 ) { // SMESH_subMesh with ID 0 may exist, don't use it!
588 rootSM.push_back( mesh.GetSubMesh( shape ));
591 for ( TopoDS_Iterator it( shape ); it.More(); it.Next() )
592 rootSM.push_back( mesh.GetSubMesh( it.Value() ));
597 // add subshapes of empty submeshes
598 list< SMESH_subMesh* >::iterator rootIt = rootSM.begin(), rootEnd = rootSM.end();
599 for ( ; rootIt != rootEnd; ++rootIt ) {
600 SMESH_subMesh * root = *rootIt;
601 SMESH_subMeshIteratorPtr smIt = root->getDependsOnIterator(/*includeSelf=*/true,
602 /*complexShapeFirst=*/true);
603 // to find a right orientation of subshapes (PAL20462)
604 TopTools_IndexedMapOfShape subShapes;
605 TopExp::MapShapes(root->GetSubShape(), subShapes);
606 while ( smIt->more() )
608 SMESH_subMesh* sm = smIt->next();
609 TopoDS_Shape shape = sm->GetSubShape();
610 totNbFaces += ( shape.ShapeType() == TopAbs_FACE );
611 if ( intern && intern->isShapeToPrecompute( shape ))
613 if ( !meshedSM || sm->IsEmpty() )
615 if ( shape.ShapeType() != TopAbs_VERTEX )
616 shape = subShapes( subShapes.FindIndex( shape ));// shape -> index -> oriented shape
617 if ( shape.Orientation() >= TopAbs_INTERNAL )
618 shape.Orientation( TopAbs_FORWARD ); // isuue 0020676
619 switch ( shape.ShapeType() ) {
620 case TopAbs_FACE : occgeo.fmap.Add( shape ); break;
621 case TopAbs_EDGE : occgeo.emap.Add( shape ); break;
622 case TopAbs_VERTEX: occgeo.vmap.Add( shape ); break;
623 case TopAbs_SOLID :occgeo.somap.Add( shape ); break;
627 // collect submeshes of meshed shapes
630 const int dim = SMESH_Gen::GetShapeDim( shape );
631 meshedSM[ dim ].push_back( sm );
635 occgeo.facemeshstatus.SetSize (totNbFaces);
636 occgeo.facemeshstatus = 0;
637 occgeo.face_maxh_modified.SetSize(totNbFaces);
638 occgeo.face_maxh_modified = 0;
639 occgeo.face_maxh.SetSize(totNbFaces);
640 occgeo.face_maxh = netgen::mparam.maxh;
643 //================================================================================
645 * \brief Return a default min size value suitable for the given geometry.
647 //================================================================================
649 double NETGENPlugin_Mesher::GetDefaultMinSize(const TopoDS_Shape& geom,
650 const double maxSize)
652 updateTriangulation( geom );
656 const int* pi[4] = { &i1, &i2, &i3, &i1 };
659 TopExp_Explorer fExp( geom, TopAbs_FACE );
660 for ( ; fExp.More(); fExp.Next() )
662 Handle(Poly_Triangulation) triangulation =
663 BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
664 if ( triangulation.IsNull() ) continue;
665 const double fTol = BRep_Tool::Tolerance( TopoDS::Face( fExp.Current() ));
666 const TColgp_Array1OfPnt& points = triangulation->Nodes();
667 const Poly_Array1OfTriangle& trias = triangulation->Triangles();
668 for ( int iT = trias.Lower(); iT <= trias.Upper(); ++iT )
670 trias(iT).Get( i1, i2, i3 );
671 for ( int j = 0; j < 3; ++j )
673 double dist2 = points(*pi[j]).SquareDistance( points( *pi[j+1] ));
674 if ( dist2 < minh && fTol*fTol < dist2 )
676 bb.Add( points(*pi[j]));
680 if ( minh > 0.25 * bb.SquareExtent() ) // simple geometry, rough triangulation
682 minh = 1e-3 * sqrt( bb.SquareExtent());
683 //cout << "BND BOX minh = " <<minh << endl;
687 minh = 3 * sqrt( minh ); // triangulation for visualization is rather fine
688 //cout << "TRIANGULATION minh = " <<minh << endl;
690 if ( minh > 0.5 * maxSize )
696 //================================================================================
698 * \brief Restrict size of elements at a given point
700 //================================================================================
702 void NETGENPlugin_Mesher::RestrictLocalSize(netgen::Mesh& ngMesh,
705 const bool overrideMinH)
707 if ( size <= std::numeric_limits<double>::min() )
709 if ( netgen::mparam.minh > size )
713 ngMesh.SetMinimalH( size );
714 netgen::mparam.minh = size;
718 size = netgen::mparam.minh;
721 netgen::Point3d pi(p.X(), p.Y(), p.Z());
722 ngMesh.RestrictLocalH( pi, size );
725 //================================================================================
727 * \brief fill ngMesh with nodes and elements of computed submeshes
729 //================================================================================
731 bool NETGENPlugin_Mesher::FillNgMesh(netgen::OCCGeometry& occgeom,
732 netgen::Mesh& ngMesh,
733 vector<const SMDS_MeshNode*>& nodeVec,
734 const list< SMESH_subMesh* > & meshedSM,
735 SMESH_MesherHelper* quadHelper,
736 SMESH_ProxyMesh::Ptr proxyMesh)
738 TNode2IdMap nodeNgIdMap;
739 for ( size_t i = 1; i < nodeVec.size(); ++i )
740 nodeNgIdMap.insert( make_pair( nodeVec[i], i ));
742 TopTools_MapOfShape visitedShapes;
743 map< SMESH_subMesh*, set< int > > visitedEdgeSM2Faces;
744 set< SMESH_subMesh* > computedSM( meshedSM.begin(), meshedSM.end() );
746 SMESH_MesherHelper helper (*_mesh);
748 int faceNgID = ngMesh.GetNFD();
750 list< SMESH_subMesh* >::const_iterator smIt, smEnd = meshedSM.end();
751 for ( smIt = meshedSM.begin(); smIt != smEnd; ++smIt )
753 SMESH_subMesh* sm = *smIt;
754 if ( !visitedShapes.Add( sm->GetSubShape() ))
757 const SMESHDS_SubMesh * smDS = sm->GetSubMeshDS();
758 if ( !smDS ) continue;
760 switch ( sm->GetSubShape().ShapeType() )
762 case TopAbs_EDGE: { // EDGE
763 // ----------------------
764 TopoDS_Edge geomEdge = TopoDS::Edge( sm->GetSubShape() );
765 if ( geomEdge.Orientation() >= TopAbs_INTERNAL )
766 geomEdge.Orientation( TopAbs_FORWARD ); // issue 0020676
768 // Add ng segments for each not meshed FACE the EDGE bounds
769 PShapeIteratorPtr fIt = helper.GetAncestors( geomEdge, *sm->GetFather(), TopAbs_FACE );
770 while ( const TopoDS_Shape * anc = fIt->next() )
772 faceNgID = occgeom.fmap.FindIndex( *anc );
774 continue; // meshed face
776 int faceSMDSId = helper.GetMeshDS()->ShapeToIndex( *anc );
777 if ( visitedEdgeSM2Faces[ sm ].count( faceSMDSId ))
778 continue; // already treated EDGE
780 TopoDS_Face face = TopoDS::Face( occgeom.fmap( faceNgID ));
781 if ( face.Orientation() >= TopAbs_INTERNAL )
782 face.Orientation( TopAbs_FORWARD ); // issue 0020676
784 // get all meshed EDGEs of the FACE connected to geomEdge (issue 0021140)
785 helper.SetSubShape( face );
786 list< TopoDS_Edge > edges = getConnectedEdges( geomEdge, face, computedSM, helper,
787 visitedEdgeSM2Faces );
789 continue; // wrong ancestor?
791 // find out orientation of <edges> within <face>
792 TopoDS_Edge eNotSeam = edges.front();
793 if ( helper.HasSeam() )
795 list< TopoDS_Edge >::iterator eIt = edges.begin();
796 while ( helper.IsRealSeam( *eIt )) ++eIt;
797 if ( eIt != edges.end() )
800 TopAbs_Orientation fOri = helper.GetSubShapeOri( face, eNotSeam );
801 bool isForwad = ( fOri == eNotSeam.Orientation() || fOri >= TopAbs_INTERNAL );
803 // get all nodes from connected <edges>
804 const bool isQuad = smDS->IsQuadratic();
805 StdMeshers_FaceSide fSide( face, edges, _mesh, isForwad, isQuad );
806 const vector<UVPtStruct>& points = fSide.GetUVPtStruct();
807 if ( points.empty() )
808 return false; // invalid node params?
809 int i, nbSeg = fSide.NbSegments();
811 // remember EDGEs of fSide to treat only once
812 for ( int iE = 0; iE < fSide.NbEdges(); ++iE )
813 visitedEdgeSM2Faces[ helper.GetMesh()->GetSubMesh( fSide.Edge(iE )) ].insert(faceSMDSId);
815 double otherSeamParam = 0;
820 int prevNgId = ngNodeId( points[0].node, ngMesh, nodeNgIdMap );
822 for ( i = 0; i < nbSeg; ++i )
824 const UVPtStruct& p1 = points[ i ];
825 const UVPtStruct& p2 = points[ i+1 ];
827 if ( p1.node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX ) //an EDGE begins
830 if ( helper.IsRealSeam( p1.node->getshapeId() ))
832 TopoDS_Edge e = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
833 isSeam = helper.IsRealSeam( e );
836 otherSeamParam = helper.GetOtherParam( helper.GetPeriodicIndex() & 1 ? p2.u : p2.v );
843 seg[1] = prevNgId = ngNodeId( p2.node, ngMesh, nodeNgIdMap );
844 // node param on curve
845 seg.epgeominfo[ 0 ].dist = p1.param;
846 seg.epgeominfo[ 1 ].dist = p2.param;
848 seg.epgeominfo[ 0 ].u = p1.u;
849 seg.epgeominfo[ 0 ].v = p1.v;
850 seg.epgeominfo[ 1 ].u = p2.u;
851 seg.epgeominfo[ 1 ].v = p2.v;
853 //geomEdge = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
854 //seg.epgeominfo[ 0 ].edgenr = seg.epgeominfo[ 1 ].edgenr = occgeom.emap.FindIndex( geomEdge );
856 //seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
857 seg.si = faceNgID; // = geom.fmap.FindIndex (face);
858 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
859 ngMesh.AddSegment (seg);
861 SMESH_TNodeXYZ np1( p1.node ), np2( p2.node );
862 RestrictLocalSize( ngMesh, 0.5*(np1+np2), (np1-np2).Modulus() );
865 cout << "Segment: " << seg.edgenr << " on SMESH face " << helper.GetMeshDS()->ShapeToIndex( face ) << endl
866 << "\tface index: " << seg.si << endl
867 << "\tp1: " << seg[0] << endl
868 << "\tp2: " << seg[1] << endl
869 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
870 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
871 //<< "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
872 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
873 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl;
874 //<< "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
878 if ( helper.GetPeriodicIndex() && 1 ) {
879 seg.epgeominfo[ 0 ].u = otherSeamParam;
880 seg.epgeominfo[ 1 ].u = otherSeamParam;
881 swap (seg.epgeominfo[0].v, seg.epgeominfo[1].v);
883 seg.epgeominfo[ 0 ].v = otherSeamParam;
884 seg.epgeominfo[ 1 ].v = otherSeamParam;
885 swap (seg.epgeominfo[0].u, seg.epgeominfo[1].u);
887 swap( seg[0], seg[1] );
888 swap( seg.epgeominfo[0].dist, seg.epgeominfo[1].dist );
889 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
890 ngMesh.AddSegment( seg );
892 cout << "Segment: " << seg.edgenr << endl
893 << "\t is SEAM (reverse) of the previous. "
894 << " Other " << (helper.GetPeriodicIndex() && 1 ? "U" : "V")
895 << " = " << otherSeamParam << endl;
898 else if ( fOri == TopAbs_INTERNAL )
900 swap( seg[0], seg[1] );
901 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
902 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
903 ngMesh.AddSegment( seg );
905 cout << "Segment: " << seg.edgenr << endl << "\t is REVERSE of the previous" << endl;
909 } // loop on geomEdge ancestors
911 if ( quadHelper ) // remember medium nodes of sub-meshes
913 SMDS_ElemIteratorPtr edges = smDS->GetElements();
914 while ( edges->more() )
916 const SMDS_MeshElement* e = edges->next();
917 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshEdge*>( e )))
923 } // case TopAbs_EDGE
925 case TopAbs_FACE: { // FACE
926 // ----------------------
927 const TopoDS_Face& geomFace = TopoDS::Face( sm->GetSubShape() );
928 helper.SetSubShape( geomFace );
929 bool isInternalFace = ( geomFace.Orientation() == TopAbs_INTERNAL );
931 // Find solids the geomFace bounds
932 int solidID1 = 0, solidID2 = 0;
933 StdMeshers_QuadToTriaAdaptor* quadAdaptor =
934 dynamic_cast<StdMeshers_QuadToTriaAdaptor*>( proxyMesh.get() );
937 solidID1 = occgeom.somap.FindIndex( quadAdaptor->GetShape() );
941 PShapeIteratorPtr solidIt = helper.GetAncestors( geomFace, *sm->GetFather(), TopAbs_SOLID);
942 while ( const TopoDS_Shape * solid = solidIt->next() )
944 int id = occgeom.somap.FindIndex ( *solid );
945 if ( solidID1 && id != solidID1 ) solidID2 = id;
949 // Add ng face descriptors of meshed faces
951 ngMesh.AddFaceDescriptor( netgen::FaceDescriptor( faceNgID, solidID1, solidID2, 0 ));
953 // if second oreder is required, even already meshed faces must be passed to NETGEN
954 int fID = occgeom.fmap.Add( geomFace );
955 if ( occgeom.facemeshstatus.Size() < fID ) occgeom.facemeshstatus.SetSize( fID );
956 occgeom.facemeshstatus[ fID-1 ] = netgen::FACE_MESHED_OK;
957 while ( fID < faceNgID ) // geomFace is already in occgeom.fmap, add a copy
959 fID = occgeom.fmap.Add( BRepBuilderAPI_Copy( geomFace, /*copyGeom=*/false ));
960 if ( occgeom.facemeshstatus.Size() < fID ) occgeom.facemeshstatus.SetSize( fID );
961 occgeom.facemeshstatus[ fID-1 ] = netgen::FACE_MESHED_OK;
963 // Problem with the second order in a quadrangular mesh remains.
964 // 1) All quadrangles generated by NETGEN are moved to an inexistent face
965 // by FillSMesh() (find "AddFaceDescriptor")
966 // 2) Temporary triangles generated by StdMeshers_QuadToTriaAdaptor
967 // are on faces where quadrangles were.
968 // Due to these 2 points, wrong geom faces are used while conversion to quadratic
969 // of the mentioned above quadrangles and triangles
971 // Orient the face correctly in solidID1 (issue 0020206)
972 bool reverse = false;
974 TopoDS_Shape solid = occgeom.somap( solidID1 );
975 TopAbs_Orientation faceOriInSolid = helper.GetSubShapeOri( solid, geomFace );
976 if ( faceOriInSolid >= 0 )
978 helper.IsReversedSubMesh( TopoDS::Face( geomFace.Oriented( faceOriInSolid )));
981 // Add surface elements
983 netgen::Element2d tri(3);
984 tri.SetIndex( faceNgID );
985 SMESH_TNodeXYZ xyz[3];
987 #ifdef DUMP_TRIANGLES
988 cout << "SMESH face " << helper.GetMeshDS()->ShapeToIndex( geomFace )
989 << " internal="<<isInternalFace << endl;
992 smDS = proxyMesh->GetSubMesh( geomFace );
994 SMDS_ElemIteratorPtr faces = smDS->GetElements();
995 while ( faces->more() )
997 const SMDS_MeshElement* f = faces->next();
998 if ( f->NbNodes() % 3 != 0 ) // not triangle
1000 PShapeIteratorPtr solidIt=helper.GetAncestors(geomFace,*sm->GetFather(),TopAbs_SOLID);
1001 if ( const TopoDS_Shape * solid = solidIt->next() )
1002 sm = _mesh->GetSubMesh( *solid );
1003 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1004 smError.reset( new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH,"Not triangle sub-mesh"));
1005 smError->myBadElements.push_back( f );
1009 for ( int i = 0; i < 3; ++i )
1011 const SMDS_MeshNode* node = f->GetNode( i ), * inFaceNode=0;
1014 // get node UV on face
1015 int shapeID = node->getshapeId();
1016 if ( helper.IsSeamShape( shapeID ))
1018 if ( helper.IsSeamShape( f->GetNodeWrap( i+1 )->getshapeId() ))
1019 inFaceNode = f->GetNodeWrap( i-1 );
1021 inFaceNode = f->GetNodeWrap( i+1 );
1023 gp_XY uv = helper.GetNodeUV( geomFace, node, inFaceNode );
1025 int ind = reverse ? 3-i : i+1;
1026 tri.GeomInfoPi(ind).u = uv.X();
1027 tri.GeomInfoPi(ind).v = uv.Y();
1028 tri.PNum (ind) = ngNodeId( node, ngMesh, nodeNgIdMap );
1031 // pass a triangle size to NG size-map
1032 double size = ( ( xyz[0] - xyz[1] ).Modulus() +
1033 ( xyz[1] - xyz[2] ).Modulus() +
1034 ( xyz[2] - xyz[0] ).Modulus() ) / 3;
1035 gp_XYZ gc = ( xyz[0] + xyz[1] + xyz[2] ) / 3;
1036 RestrictLocalSize( ngMesh, gc, size, /*overrideMinH=*/false );
1038 ngMesh.AddSurfaceElement (tri);
1039 #ifdef DUMP_TRIANGLES
1040 cout << tri << endl;
1043 if ( isInternalFace )
1045 swap( tri[1], tri[2] );
1046 ngMesh.AddSurfaceElement (tri);
1047 #ifdef DUMP_TRIANGLES
1048 cout << tri << endl;
1053 if ( quadHelper ) // remember medium nodes of sub-meshes
1055 SMDS_ElemIteratorPtr faces = smDS->GetElements();
1056 while ( faces->more() )
1058 const SMDS_MeshElement* f = faces->next();
1059 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshFace*>( f )))
1065 } // case TopAbs_FACE
1067 case TopAbs_VERTEX: { // VERTEX
1068 // --------------------------
1069 // issue 0021405. Add node only if a VERTEX is shared by a not meshed EDGE,
1070 // else netgen removes a free node and nodeVector becomes invalid
1071 PShapeIteratorPtr ansIt = helper.GetAncestors( sm->GetSubShape(),
1075 while ( const TopoDS_Shape* e = ansIt->next() )
1077 SMESH_subMesh* eSub = helper.GetMesh()->GetSubMesh( *e );
1078 if (( toAdd = ( eSub->IsEmpty() && !SMESH_Algo::isDegenerated( TopoDS::Edge( *e )))))
1083 SMDS_NodeIteratorPtr nodeIt = smDS->GetNodes();
1084 if ( nodeIt->more() )
1085 ngNodeId( nodeIt->next(), ngMesh, nodeNgIdMap );
1091 } // loop on submeshes
1094 nodeVec.resize( ngMesh.GetNP() + 1 );
1095 TNode2IdMap::iterator node_NgId, nodeNgIdEnd = nodeNgIdMap.end();
1096 for ( node_NgId = nodeNgIdMap.begin(); node_NgId != nodeNgIdEnd; ++node_NgId)
1097 nodeVec[ node_NgId->second ] = node_NgId->first;
1102 //================================================================================
1104 * \brief Duplicate mesh faces on internal geom faces
1106 //================================================================================
1108 void NETGENPlugin_Mesher::FixIntFaces(const netgen::OCCGeometry& occgeom,
1109 netgen::Mesh& ngMesh,
1110 NETGENPlugin_Internals& internalShapes)
1112 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1114 // find ng indices of internal faces
1116 for ( int ngFaceID = 1; ngFaceID <= occgeom.fmap.Extent(); ++ngFaceID )
1118 int smeshID = meshDS->ShapeToIndex( occgeom.fmap( ngFaceID ));
1119 if ( internalShapes.isInternalShape( smeshID ))
1120 ngFaceIds.insert( ngFaceID );
1122 if ( !ngFaceIds.empty() )
1125 int i, nbFaces = ngMesh.GetNSE();
1126 for ( i = 1; i <= nbFaces; ++i)
1128 netgen::Element2d elem = ngMesh.SurfaceElement(i);
1129 if ( ngFaceIds.count( elem.GetIndex() ))
1131 swap( elem[1], elem[2] );
1132 ngMesh.AddSurfaceElement (elem);
1138 //================================================================================
1140 * \brief Tries to heal the mesh on a FACE. The FACE is supposed to be partially
1141 * meshed due to NETGEN failure
1142 * \param [in] occgeom - geometry
1143 * \param [in,out] ngMesh - the mesh to fix
1144 * \param [inout] faceID - ID of the FACE to fix the mesh on
1145 * \return bool - is mesh is or becomes OK
1147 //================================================================================
1149 bool NETGENPlugin_Mesher::FixFaceMesh(const netgen::OCCGeometry& occgeom,
1150 netgen::Mesh& ngMesh,
1153 // we address a case where the FACE is almost fully meshed except small holes
1154 // of usually triangular shape at FACE boundary (IPAL52861)
1156 // The case appeared to be not simple: holes only look triangular but
1157 // indeed are a self intersecting polygon. A reason of the bug was in coincident
1158 // NG points on a seam edge. But the code below is very nice, leave it for
1163 if ( occgeom.fmap.Extent() < faceID )
1165 //const TopoDS_Face& face = TopoDS::Face( occgeom.fmap( faceID ));
1167 // find free links on the FACE
1168 NCollection_Map<Link> linkMap;
1169 for ( int iF = 1; iF <= ngMesh.GetNSE(); ++iF )
1171 const netgen::Element2d& elem = ngMesh.SurfaceElement(iF);
1172 if ( faceID != elem.GetIndex() )
1174 int n0 = elem[ elem.GetNP() - 1 ];
1175 for ( int i = 0; i < elem.GetNP(); ++i )
1178 Link link( n0, n1 );
1179 if ( !linkMap.Add( link ))
1180 linkMap.Remove( link );
1184 // add/remove boundary links
1185 for ( int iSeg = 1; iSeg <= ngMesh.GetNSeg(); ++iSeg )
1187 const netgen::Segment& seg = ngMesh.LineSegment( iSeg );
1188 if ( seg.si != faceID ) // !edgeIDs.Contains( seg.edgenr ))
1190 Link link( seg[1], seg[0] ); // reverse!!!
1191 if ( !linkMap.Add( link ))
1192 linkMap.Remove( link );
1194 if ( linkMap.IsEmpty() )
1196 if ( linkMap.Extent() < 3 )
1199 // make triangles of the links
1201 netgen::Element2d tri(3);
1202 tri.SetIndex ( faceID );
1204 NCollection_Map<Link>::Iterator linkIt( linkMap );
1205 Link link1 = linkIt.Value();
1206 // look for a link connected to link1
1207 NCollection_Map<Link>::Iterator linkIt2 = linkIt;
1208 for ( linkIt2.Next(); linkIt2.More(); linkIt2.Next() )
1210 const Link& link2 = linkIt2.Value();
1211 if ( link2.IsConnected( link1 ))
1213 // look for a link connected to both link1 and link2
1214 NCollection_Map<Link>::Iterator linkIt3 = linkIt2;
1215 for ( linkIt3.Next(); linkIt3.More(); linkIt3.Next() )
1217 const Link& link3 = linkIt3.Value();
1218 if ( link3.IsConnected( link1 ) &&
1219 link3.IsConnected( link2 ) )
1224 tri[2] = ( link2.Contains( link1.n1 ) ? link2.n1 : link3.n1 );
1225 if ( tri[0] == tri[2] || tri[1] == tri[2] )
1227 ngMesh.AddSurfaceElement( tri );
1229 // prepare for the next tria search
1230 if ( linkMap.Extent() == 3 )
1232 linkMap.Remove( link3 );
1233 linkMap.Remove( link2 );
1235 linkMap.Remove( link1 );
1236 link1 = linkIt.Value();
1249 //================================================================================
1250 // define gp_XY_Subtracted pointer to function calling gp_XY::Subtracted(gp_XY)
1251 gp_XY_FunPtr(Subtracted);
1252 //gp_XY_FunPtr(Added);
1254 //================================================================================
1256 * \brief Evaluate distance between two 2d points along the surface
1258 //================================================================================
1260 double evalDist( const gp_XY& uv1,
1262 const Handle(Geom_Surface)& surf,
1263 const int stopHandler=-1)
1265 if ( stopHandler > 0 ) // continue recursion
1267 gp_XY mid = SMESH_MesherHelper::GetMiddleUV( surf, uv1, uv2 );
1268 return evalDist( uv1,mid, surf, stopHandler-1 ) + evalDist( mid,uv2, surf, stopHandler-1 );
1270 double dist3D = surf->Value( uv1.X(), uv1.Y() ).Distance( surf->Value( uv2.X(), uv2.Y() ));
1271 if ( stopHandler == 0 ) // stop recursion
1274 // start recursion if necessary
1275 double dist2D = SMESH_MesherHelper::ApplyIn2D(surf, uv1, uv2, gp_XY_Subtracted, 0).Modulus();
1276 if ( fabs( dist3D - dist2D ) < dist2D * 1e-10 )
1277 return dist3D; // equal parametrization of a planar surface
1279 return evalDist( uv1, uv2, surf, 3 ); // start recursion
1282 //================================================================================
1284 * \brief Data of vertex internal in geom face
1286 //================================================================================
1290 gp_XY uv; //!< UV in face parametric space
1291 int ngId; //!< ng id of corrsponding node
1292 gp_XY uvClose; //!< UV of closest boundary node
1293 int ngIdClose; //!< ng id of closest boundary node
1296 //================================================================================
1298 * \brief Data of vertex internal in solid
1300 //================================================================================
1304 int ngId; //!< ng id of corresponding node
1305 int ngIdClose; //!< ng id of closest 2d mesh element
1306 int ngIdCloseN; //!< ng id of closest node of the closest 2d mesh element
1309 inline double dist2(const netgen::MeshPoint& p1, const netgen::MeshPoint& p2)
1311 return gp_Pnt( NGPOINT_COORDS(p1)).SquareDistance( gp_Pnt( NGPOINT_COORDS(p2)));
1315 //================================================================================
1317 * \brief Make netgen take internal vertices in faces into account by adding
1318 * segments including internal vertices
1320 * This function works in supposition that 1D mesh is already computed in ngMesh
1322 //================================================================================
1324 void NETGENPlugin_Mesher::AddIntVerticesInFaces(const netgen::OCCGeometry& occgeom,
1325 netgen::Mesh& ngMesh,
1326 vector<const SMDS_MeshNode*>& nodeVec,
1327 NETGENPlugin_Internals& internalShapes)
1329 if ((int) nodeVec.size() < ngMesh.GetNP() )
1330 nodeVec.resize( ngMesh.GetNP(), 0 );
1332 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1333 SMESH_MesherHelper helper( internalShapes.getMesh() );
1335 const map<int,list<int> >& face2Vert = internalShapes.getFacesWithVertices();
1336 map<int,list<int> >::const_iterator f2v = face2Vert.begin();
1337 for ( ; f2v != face2Vert.end(); ++f2v )
1339 const TopoDS_Face& face = TopoDS::Face( meshDS->IndexToShape( f2v->first ));
1340 if ( face.IsNull() ) continue;
1341 int faceNgID = occgeom.fmap.FindIndex (face);
1342 if ( faceNgID < 0 ) continue;
1344 TopLoc_Location loc;
1345 Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc);
1347 helper.SetSubShape( face );
1348 helper.SetElementsOnShape( true );
1350 // Get data of internal vertices and add them to ngMesh
1352 multimap< double, TIntVData > dist2VData; // sort vertices by distance from boundary nodes
1354 int i, nbSegInit = ngMesh.GetNSeg();
1356 // boundary characteristics
1357 double totSegLen2D = 0;
1360 const list<int>& iVertices = f2v->second;
1361 list<int>::const_iterator iv = iVertices.begin();
1362 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1365 // get node on vertex
1366 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1367 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1370 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1371 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1372 nV = SMESH_Algo::VertexNode( V, meshDS );
1373 if ( !nV ) continue;
1376 netgen::MeshPoint mp( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1377 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1378 vData.ngId = ngMesh.GetNP();
1379 nodeVec.push_back( nV );
1383 vData.uv = helper.GetNodeUV( face, nV, 0, &uvOK );
1384 if ( !uvOK ) helper.CheckNodeUV( face, nV, vData.uv, BRep_Tool::Tolerance(V),/*force=*/1);
1386 // loop on all segments of the face to find the node closest to vertex and to count
1387 // average segment 2d length
1388 double closeDist2 = numeric_limits<double>::max(), dist2;
1390 for (i = 1; i <= ngMesh.GetNSeg(); ++i)
1392 netgen::Segment & seg = ngMesh.LineSegment(i);
1393 if ( seg.si != faceNgID ) continue;
1395 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1397 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1398 if ( ngIdLast == seg[ iEnd ] ) continue;
1399 dist2 = helper.ApplyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1400 if ( dist2 < closeDist2 )
1401 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1402 ngIdLast = seg[ iEnd ];
1406 totSegLen2D += helper.ApplyIn2D(surf, uv[0], uv[1], gp_XY_Subtracted, false).Modulus();
1410 dist2VData.insert( make_pair( closeDist2, vData ));
1413 if ( totNbSeg == 0 ) break;
1414 double avgSegLen2d = totSegLen2D / totNbSeg;
1416 // Loop on vertices to add segments
1418 multimap< double, TIntVData >::iterator dist_vData = dist2VData.begin();
1419 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1421 double closeDist2 = dist_vData->first, dist2;
1422 TIntVData & vData = dist_vData->second;
1424 // try to find more close node among segments added for internal vertices
1425 for (i = nbSegInit+1; i <= ngMesh.GetNSeg(); ++i)
1427 netgen::Segment & seg = ngMesh.LineSegment(i);
1428 if ( seg.si != faceNgID ) continue;
1430 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1432 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1433 dist2 = helper.ApplyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1434 if ( dist2 < closeDist2 )
1435 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1438 // decide whether to use the closest node as the second end of segment or to
1439 // create a new point
1440 int segEnd1 = vData.ngId;
1441 int segEnd2 = vData.ngIdClose; // to use closest node
1442 gp_XY uvV = vData.uv, uvP = vData.uvClose;
1443 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1444 double nodeDist2D = sqrt( closeDist2 );
1445 double nodeDist3D = evalDist( vData.uv, vData.uvClose, surf );
1446 bool avgLenOK = ( avgSegLen2d < 0.75 * nodeDist2D );
1447 bool hintLenOK = ( segLenHint < 0.75 * nodeDist3D );
1448 //cout << "uvV " << uvV.X() <<","<<uvV.Y() << " ";
1449 if ( hintLenOK || avgLenOK )
1451 // create a point between the closest node and V
1454 double r = min( 0.5, ( hintLenOK ? segLenHint/nodeDist3D : avgSegLen2d/nodeDist2D ));
1455 // direction from V to closet node in 2D
1456 gp_Dir2d v2n( helper.ApplyIn2D(surf, uvP, uvV, gp_XY_Subtracted, false ));
1458 uvP = vData.uv + r * nodeDist2D * v2n.XY();
1459 gp_Pnt P = surf->Value( uvP.X(), uvP.Y() ).Transformed( loc );
1461 netgen::MeshPoint mp( netgen::Point<3> (P.X(), P.Y(), P.Z()));
1462 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1463 segEnd2 = ngMesh.GetNP();
1464 //cout << "Middle " << r << " uv " << uvP.X() << "," << uvP.Y() << "( " << ngMesh.Point(segEnd2).X()<<","<<ngMesh.Point(segEnd2).Y()<<","<<ngMesh.Point(segEnd2).Z()<<" )"<< endl;
1465 SMDS_MeshNode * nP = helper.AddNode(P.X(), P.Y(), P.Z());
1466 nodeVec.push_back( nP );
1468 //else cout << "at Node " << " uv " << uvP.X() << "," << uvP.Y() << endl;
1471 netgen::Segment seg;
1473 if ( segEnd1 > segEnd2 ) swap( segEnd1, segEnd2 ), swap( uvV, uvP );
1474 seg[0] = segEnd1; // ng node id
1475 seg[1] = segEnd2; // ng node id
1476 seg.edgenr = ngMesh.GetNSeg() + 1;// segment id
1479 seg.epgeominfo[ 0 ].dist = 0; // param on curve
1480 seg.epgeominfo[ 0 ].u = uvV.X();
1481 seg.epgeominfo[ 0 ].v = uvV.Y();
1482 seg.epgeominfo[ 1 ].dist = 1; // param on curve
1483 seg.epgeominfo[ 1 ].u = uvP.X();
1484 seg.epgeominfo[ 1 ].v = uvP.Y();
1486 // seg.epgeominfo[ 0 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1487 // seg.epgeominfo[ 1 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1489 ngMesh.AddSegment (seg);
1491 // add reverse segment
1492 swap( seg[0], seg[1] );
1493 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1494 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1495 ngMesh.AddSegment (seg);
1501 //================================================================================
1503 * \brief Make netgen take internal vertices in solids into account by adding
1504 * faces including internal vertices
1506 * This function works in supposition that 2D mesh is already computed in ngMesh
1508 //================================================================================
1510 void NETGENPlugin_Mesher::AddIntVerticesInSolids(const netgen::OCCGeometry& occgeom,
1511 netgen::Mesh& ngMesh,
1512 vector<const SMDS_MeshNode*>& nodeVec,
1513 NETGENPlugin_Internals& internalShapes)
1515 #ifdef DUMP_TRIANGLES_SCRIPT
1516 // create a python script making a mesh containing triangles added for internal vertices
1517 ofstream py(DUMP_TRIANGLES_SCRIPT);
1518 py << "import SMESH"<< endl
1519 << "from salome.smesh import smeshBuilder"<<endl
1520 << "smesh = smeshBuilder.New(salome.myStudy)"<<endl
1521 << "m = smesh.Mesh(name='triangles')" << endl;
1523 if ((int) nodeVec.size() < ngMesh.GetNP() )
1524 nodeVec.resize( ngMesh.GetNP(), 0 );
1526 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1527 SMESH_MesherHelper helper( internalShapes.getMesh() );
1529 const map<int,list<int> >& so2Vert = internalShapes.getSolidsWithVertices();
1530 map<int,list<int> >::const_iterator s2v = so2Vert.begin();
1531 for ( ; s2v != so2Vert.end(); ++s2v )
1533 const TopoDS_Shape& solid = meshDS->IndexToShape( s2v->first );
1534 if ( solid.IsNull() ) continue;
1535 int solidNgID = occgeom.somap.FindIndex (solid);
1536 if ( solidNgID < 0 && !occgeom.somap.IsEmpty() ) continue;
1538 helper.SetSubShape( solid );
1539 helper.SetElementsOnShape( true );
1541 // find ng indices of faces within the solid
1543 for (TopExp_Explorer fExp(solid, TopAbs_FACE); fExp.More(); fExp.Next() )
1544 ngFaceIds.insert( occgeom.fmap.FindIndex( fExp.Current() ));
1545 if ( ngFaceIds.size() == 1 && *ngFaceIds.begin() == 0 )
1546 ngFaceIds.insert( 1 );
1548 // Get data of internal vertices and add them to ngMesh
1550 multimap< double, TIntVSoData > dist2VData; // sort vertices by distance from ng faces
1552 int i, nbFaceInit = ngMesh.GetNSE();
1554 // boundary characteristics
1555 double totSegLen = 0;
1558 const list<int>& iVertices = s2v->second;
1559 list<int>::const_iterator iv = iVertices.begin();
1560 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1563 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1565 // get node on vertex
1566 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1569 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1570 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1571 nV = SMESH_Algo::VertexNode( V, meshDS );
1572 if ( !nV ) continue;
1575 netgen::MeshPoint mpV( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1576 ngMesh.AddPoint ( mpV, 1, netgen::FIXEDPOINT );
1577 vData.ngId = ngMesh.GetNP();
1578 nodeVec.push_back( nV );
1580 // loop on all 2d elements to find the one closest to vertex and to count
1581 // average segment length
1582 double closeDist2 = numeric_limits<double>::max(), avgDist2;
1583 for (i = 1; i <= ngMesh.GetNSE(); ++i)
1585 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1586 if ( !ngFaceIds.count( elem.GetIndex() )) continue;
1588 multimap< double, int> dist2nID; // sort nodes of element by distance from V
1589 for ( int j = 0; j < elem.GetNP(); ++j)
1591 netgen::MeshPoint mp = ngMesh.Point( elem[j] );
1592 double d2 = dist2( mpV, mp );
1593 dist2nID.insert( make_pair( d2, elem[j] ));
1594 avgDist2 += d2 / elem.GetNP();
1596 totNbSeg++, totSegLen+= sqrt( dist2( mp, ngMesh.Point( elem[(j+1)%elem.GetNP()])));
1598 double dist = dist2nID.begin()->first; //avgDist2;
1599 if ( dist < closeDist2 )
1600 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= dist;
1602 dist2VData.insert( make_pair( closeDist2, vData ));
1605 if ( totNbSeg == 0 ) break;
1606 double avgSegLen = totSegLen / totNbSeg;
1608 // Loop on vertices to add triangles
1610 multimap< double, TIntVSoData >::iterator dist_vData = dist2VData.begin();
1611 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1613 double closeDist2 = dist_vData->first;
1614 TIntVSoData & vData = dist_vData->second;
1616 const netgen::MeshPoint& mpV = ngMesh.Point( vData.ngId );
1618 // try to find more close face among ones added for internal vertices
1619 for (i = nbFaceInit+1; i <= ngMesh.GetNSE(); ++i)
1621 double avgDist2 = 0;
1622 multimap< double, int> dist2nID;
1623 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1624 for ( int j = 0; j < elem.GetNP(); ++j)
1626 double d = dist2( mpV, ngMesh.Point( elem[j] ));
1627 dist2nID.insert( make_pair( d, elem[j] ));
1628 avgDist2 += d / elem.GetNP();
1629 if ( avgDist2 < closeDist2 )
1630 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= avgDist2;
1633 // sort nodes of the closest face by angle with vector from V to the closest node
1634 const double tol = numeric_limits<double>::min();
1635 map< double, int > angle2ID;
1636 const netgen::Element2d& closeFace = ngMesh.SurfaceElement( vData.ngIdClose );
1637 netgen::MeshPoint mp[2];
1638 mp[0] = ngMesh.Point( vData.ngIdCloseN );
1639 gp_XYZ p1( NGPOINT_COORDS( mp[0] ));
1640 gp_XYZ pV( NGPOINT_COORDS( mpV ));
1641 gp_Vec v2p1( pV, p1 );
1642 double distN1 = v2p1.Magnitude();
1643 if ( distN1 <= tol ) continue;
1645 for ( int j = 0; j < closeFace.GetNP(); ++j)
1647 mp[1] = ngMesh.Point( closeFace[j] );
1648 gp_Vec v2p( pV, gp_Pnt( NGPOINT_COORDS( mp[1] )) );
1649 angle2ID.insert( make_pair( v2p1.Angle( v2p ), closeFace[j]));
1651 // get node with angle of 60 degrees or greater
1652 map< double, int >::iterator angle_id = angle2ID.lower_bound( 60. * M_PI / 180. );
1653 if ( angle_id == angle2ID.end() ) angle_id = --angle2ID.end();
1654 const double minAngle = 30. * M_PI / 180.;
1655 const double angle = angle_id->first;
1656 bool angleOK = ( angle > minAngle );
1658 // find points to create a triangle
1659 netgen::Element2d tri(3);
1661 tri[0] = vData.ngId;
1662 tri[1] = vData.ngIdCloseN; // to use the closest nodes
1663 tri[2] = angle_id->second; // to use the node with best angle
1665 // decide whether to use the closest node and the node with best angle or to create new ones
1666 for ( int isBestAngleN = 0; isBestAngleN < 2; ++isBestAngleN )
1668 bool createNew = !angleOK; //, distOK = true;
1670 int triInd = isBestAngleN ? 2 : 1;
1671 mp[isBestAngleN] = ngMesh.Point( tri[triInd] );
1676 double distN2 = sqrt( dist2( mpV, mp[isBestAngleN]));
1677 createNew = ( fabs( distN2 - distN1 ) > 0.25 * distN1 );
1679 else if ( angle < tol )
1681 v2p1.SetX( v2p1.X() + 1e-3 );
1687 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1688 bool avgLenOK = ( avgSegLen < 0.75 * distN1 );
1689 bool hintLenOK = ( segLenHint < 0.75 * distN1 );
1690 createNew = (createNew || avgLenOK || hintLenOK );
1691 // we create a new node not closer than 0.5 to the closest face
1692 // in order not to clash with other close face
1693 double r = min( 0.5, ( hintLenOK ? segLenHint : avgSegLen ) / distN1 );
1694 distFromV = r * distN1;
1698 // create a new point, between the node and the vertex if angleOK
1699 gp_XYZ p( NGPOINT_COORDS( mp[isBestAngleN] ));
1700 gp_Vec v2p( pV, p ); v2p.Normalize();
1701 if ( isBestAngleN && !angleOK )
1702 p = p1 + gp_Dir( v2p.XYZ() - v2p1.XYZ()).XYZ() * distN1 * 0.95;
1704 p = pV + v2p.XYZ() * distFromV;
1706 if ( !isBestAngleN ) p1 = p, distN1 = distFromV;
1708 mp[isBestAngleN].SetPoint( netgen::Point<3> (p.X(), p.Y(), p.Z()));
1709 ngMesh.AddPoint ( mp[isBestAngleN], 1, netgen::SURFACEPOINT );
1710 tri[triInd] = ngMesh.GetNP();
1711 nodeVec.push_back( helper.AddNode( p.X(), p.Y(), p.Z()) );
1714 ngMesh.AddSurfaceElement (tri);
1715 swap( tri[1], tri[2] );
1716 ngMesh.AddSurfaceElement (tri);
1718 #ifdef DUMP_TRIANGLES_SCRIPT
1719 py << "n1 = m.AddNode( "<< mpV(0)<<", "<< mpV(1)<<", "<< mpV(2)<<") "<< endl
1720 << "n2 = m.AddNode( "<< mp[0](0)<<", "<< mp[0](1)<<", "<< mp[0](2)<<") "<< endl
1721 << "n3 = m.AddNode( "<< mp[1](0)<<", "<< mp[1](1)<<", "<< mp[1](2)<<" )" << endl
1722 << "m.AddFace([n1,n2,n3])" << endl;
1724 } // loop on internal vertices of a solid
1726 } // loop on solids with internal vertices
1729 //================================================================================
1731 * \brief Fill netgen mesh with segments of a FACE
1732 * \param ngMesh - netgen mesh
1733 * \param geom - container of OCCT geometry to mesh
1734 * \param wires - data of nodes on FACE boundary
1735 * \param helper - mesher helper holding the FACE
1736 * \param nodeVec - vector of nodes in which node index == netgen ID
1737 * \retval SMESH_ComputeErrorPtr - error description
1739 //================================================================================
1741 SMESH_ComputeErrorPtr
1742 NETGENPlugin_Mesher::AddSegmentsToMesh(netgen::Mesh& ngMesh,
1743 netgen::OCCGeometry& geom,
1744 const TSideVector& wires,
1745 SMESH_MesherHelper& helper,
1746 vector< const SMDS_MeshNode* > & nodeVec,
1747 const bool overrideMinH)
1749 // ----------------------------
1750 // Check wires and count nodes
1751 // ----------------------------
1753 for ( size_t iW = 0; iW < wires.size(); ++iW )
1755 StdMeshers_FaceSidePtr wire = wires[ iW ];
1756 if ( wire->MissVertexNode() )
1758 // Commented for issue 0020960. It worked for the case, let's wait for case where it doesn't.
1759 // It seems that there is no reason for this limitation
1761 // (new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH, "Missing nodes on vertices"));
1763 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1764 if ((int) uvPtVec.size() != wire->NbPoints() )
1765 return SMESH_ComputeError::New(COMPERR_BAD_INPUT_MESH,
1766 SMESH_Comment("Unexpected nb of points on wire ") << iW
1767 << ": " << uvPtVec.size()<<" != "<<wire->NbPoints());
1768 nbNodes += wire->NbPoints();
1770 nodeVec.reserve( nodeVec.size() + nbNodes + 1 );
1771 if ( nodeVec.empty() )
1772 nodeVec.push_back( 0 );
1774 // -----------------
1776 // -----------------
1778 const bool wasNgMeshEmpty = ( ngMesh.GetNP() < 1 ); /* true => this method is called by
1779 NETGENPlugin_NETGEN_2D_ONLY */
1781 // map for nodes on vertices since they can be shared between wires
1782 // ( issue 0020676, face_int_box.brep) and nodes built by NETGEN
1783 map<const SMDS_MeshNode*, int > node2ngID;
1784 if ( !wasNgMeshEmpty ) // fill node2ngID with nodes built by NETGEN
1786 set< int > subIDs; // ids of sub-shapes of the FACE
1787 for ( size_t iW = 0; iW < wires.size(); ++iW )
1789 StdMeshers_FaceSidePtr wire = wires[ iW ];
1790 for ( int iE = 0, nbE = wire->NbEdges(); iE < nbE; ++iE )
1792 subIDs.insert( wire->EdgeID( iE ));
1793 subIDs.insert( helper.GetMeshDS()->ShapeToIndex( wire->FirstVertex( iE )));
1796 for ( size_t ngID = 1; ngID < nodeVec.size(); ++ngID )
1797 if ( subIDs.count( nodeVec[ngID]->getshapeId() ))
1798 node2ngID.insert( make_pair( nodeVec[ngID], ngID ));
1801 const int solidID = 0, faceID = geom.fmap.FindIndex( helper.GetSubShape() );
1802 if ( ngMesh.GetNFD() < 1 )
1803 ngMesh.AddFaceDescriptor( netgen::FaceDescriptor( faceID, solidID, solidID, 0 ));
1805 for ( size_t iW = 0; iW < wires.size(); ++iW )
1807 StdMeshers_FaceSidePtr wire = wires[ iW ];
1808 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1809 const int nbSegments = wire->NbPoints() - 1;
1811 // assure the 1st node to be in node2ngID, which is needed to correctly
1812 // "close chain of segments" (see below) in case if the 1st node is not
1813 // onVertex because it is on a Viscous layer
1814 node2ngID.insert( make_pair( uvPtVec[ 0 ].node, ngMesh.GetNP() + 1 ));
1816 // compute length of every segment
1817 vector<double> segLen( nbSegments );
1818 for ( int i = 0; i < nbSegments; ++i )
1819 segLen[i] = SMESH_TNodeXYZ( uvPtVec[ i ].node ).Distance( uvPtVec[ i+1 ].node );
1821 int edgeID = 1, posID = -2;
1822 bool isInternalWire = false;
1823 double vertexNormPar = 0;
1824 const int prevNbNGSeg = ngMesh.GetNSeg();
1825 for ( int i = 0; i < nbSegments; ++i ) // loop on segments
1827 // Add the first point of a segment
1829 const SMDS_MeshNode * n = uvPtVec[ i ].node;
1830 const int posShapeID = n->getshapeId();
1831 bool onVertex = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX );
1832 bool onEdge = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_EDGE );
1834 // skip nodes on degenerated edges
1835 if ( helper.IsDegenShape( posShapeID ) &&
1836 helper.IsDegenShape( uvPtVec[ i+1 ].node->getshapeId() ))
1839 int ngID1 = ngMesh.GetNP() + 1, ngID2 = ngID1+1;
1840 if ( onVertex || ( !wasNgMeshEmpty && onEdge ) || helper.IsRealSeam( posShapeID ))
1841 ngID1 = node2ngID.insert( make_pair( n, ngID1 )).first->second;
1842 if ( ngID1 > ngMesh.GetNP() )
1844 netgen::MeshPoint mp( netgen::Point<3> (n->X(), n->Y(), n->Z()) );
1845 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1846 nodeVec.push_back( n );
1848 else // n is in ngMesh already, and ngID2 in prev segment is wrong
1850 ngID2 = ngMesh.GetNP() + 1;
1851 if ( i > 0 ) // prev segment belongs to same wire
1853 netgen::Segment& prevSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
1860 netgen::Segment seg;
1862 seg[0] = ngID1; // ng node id
1863 seg[1] = ngID2; // ng node id
1864 seg.edgenr = ngMesh.GetNSeg() + 1; // ng segment id
1865 seg.si = faceID; // = geom.fmap.FindIndex (face);
1867 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1869 const UVPtStruct& pnt = uvPtVec[ i + iEnd ];
1871 seg.epgeominfo[ iEnd ].dist = pnt.param; // param on curve
1872 seg.epgeominfo[ iEnd ].u = pnt.u;
1873 seg.epgeominfo[ iEnd ].v = pnt.v;
1875 // find out edge id and node parameter on edge
1876 onVertex = ( pnt.normParam + 1e-10 > vertexNormPar );
1877 if ( onVertex || posShapeID != posID )
1880 double normParam = pnt.normParam;
1882 normParam = 0.5 * ( uvPtVec[ i ].normParam + uvPtVec[ i+1 ].normParam );
1883 int edgeIndexInWire = wire->EdgeIndex( normParam );
1884 vertexNormPar = wire->LastParameter( edgeIndexInWire );
1885 const TopoDS_Edge& edge = wire->Edge( edgeIndexInWire );
1886 edgeID = geom.emap.FindIndex( edge );
1888 isInternalWire = ( edge.Orientation() == TopAbs_INTERNAL );
1889 // if ( onVertex ) // param on curve is different on each of two edges
1890 // seg.epgeominfo[ iEnd ].dist = helper.GetNodeU( edge, pnt.node );
1892 seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
1895 ngMesh.AddSegment (seg);
1897 // restrict size of elements near the segment
1898 SMESH_TNodeXYZ np1( n ), np2( uvPtVec[ i+1 ].node );
1899 // get an average size of adjacent segments to avoid sharp change of
1900 // element size (regression on issue 0020452, note 0010898)
1901 int iPrev = SMESH_MesherHelper::WrapIndex( i-1, nbSegments );
1902 int iNext = SMESH_MesherHelper::WrapIndex( i+1, nbSegments );
1903 double sumH = segLen[ iPrev ] + segLen[ i ] + segLen[ iNext ];
1904 int nbSeg = ( int( segLen[ iPrev ] > sumH / 100.) +
1905 int( segLen[ i ] > sumH / 100.) +
1906 int( segLen[ iNext ] > sumH / 100.));
1908 RestrictLocalSize( ngMesh, 0.5*(np1+np2), sumH / nbSeg, overrideMinH );
1910 if ( isInternalWire )
1912 swap (seg[0], seg[1]);
1913 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1914 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1915 ngMesh.AddSegment (seg);
1917 } // loop on segments on a wire
1919 // close chain of segments
1920 if ( nbSegments > 0 )
1922 netgen::Segment& lastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() - int( isInternalWire ));
1923 const SMDS_MeshNode * lastNode = uvPtVec.back().node;
1924 lastSeg[1] = node2ngID.insert( make_pair( lastNode, lastSeg[1] )).first->second;
1925 if ( lastSeg[1] > ngMesh.GetNP() )
1927 netgen::MeshPoint mp( netgen::Point<3> (lastNode->X(), lastNode->Y(), lastNode->Z()) );
1928 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1929 nodeVec.push_back( lastNode );
1931 if ( isInternalWire )
1933 netgen::Segment& realLastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
1934 realLastSeg[0] = lastSeg[1];
1938 #ifdef DUMP_SEGMENTS
1939 cout << "BEGIN WIRE " << iW << endl;
1940 for ( int i = prevNbNGSeg+1; i <= ngMesh.GetNSeg(); ++i )
1942 netgen::Segment& seg = ngMesh.LineSegment( i );
1944 netgen::Segment& prevSeg = ngMesh.LineSegment( i-1 );
1945 if ( seg[0] == prevSeg[1] && seg[1] == prevSeg[0] )
1947 cout << "Segment: " << seg.edgenr << endl << "\tis REVRESE of the previous one" << endl;
1951 cout << "Segment: " << seg.edgenr << endl
1952 << "\tp1: " << seg[0] << " n" << nodeVec[ seg[0]]->GetID() << endl
1953 << "\tp2: " << seg[1] << " n" << nodeVec[ seg[1]]->GetID() << endl
1954 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
1955 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
1956 << "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
1957 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
1958 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl
1959 << "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
1961 cout << "--END WIRE " << iW << endl;
1963 SMESH_Comment __not_unused_variable( prevNbNGSeg );
1966 } // loop on WIREs of a FACE
1968 // add a segment instead of an internal vertex
1969 if ( wasNgMeshEmpty )
1971 NETGENPlugin_Internals intShapes( *helper.GetMesh(), helper.GetSubShape(), /*is3D=*/false );
1972 AddIntVerticesInFaces( geom, ngMesh, nodeVec, intShapes );
1974 ngMesh.CalcSurfacesOfNode();
1979 //================================================================================
1981 * \brief Fill SMESH mesh according to contents of netgen mesh
1982 * \param occgeo - container of OCCT geometry to mesh
1983 * \param ngMesh - netgen mesh
1984 * \param initState - bn of entities in netgen mesh before computing
1985 * \param sMesh - SMESH mesh to fill in
1986 * \param nodeVec - vector of nodes in which node index == netgen ID
1987 * \param comment - returns problem description
1988 * \param quadHelper - holder of medium nodes of sub-meshes
1989 * \retval int - error
1991 //================================================================================
1993 int NETGENPlugin_Mesher::FillSMesh(const netgen::OCCGeometry& occgeo,
1994 netgen::Mesh& ngMesh,
1995 const NETGENPlugin_ngMeshInfo& initState,
1997 std::vector<const SMDS_MeshNode*>& nodeVec,
1998 SMESH_Comment& comment,
1999 SMESH_MesherHelper* quadHelper)
2001 int nbNod = ngMesh.GetNP();
2002 int nbSeg = ngMesh.GetNSeg();
2003 int nbFac = ngMesh.GetNSE();
2004 int nbVol = ngMesh.GetNE();
2006 SMESHDS_Mesh* meshDS = sMesh.GetMeshDS();
2008 // quadHelper is used for either
2009 // 1) making quadratic elements when a lower dimention mesh is loaded
2010 // to SMESH before convertion to quadratic by NETGEN
2011 // 2) sewing of quadratic elements with quadratic elements of sub-meshes
2012 if ( quadHelper && !quadHelper->GetIsQuadratic() && quadHelper->GetTLinkNodeMap().empty() )
2015 // -------------------------------------
2016 // Create and insert nodes into nodeVec
2017 // -------------------------------------
2019 nodeVec.resize( nbNod + 1 );
2020 int i, nbInitNod = initState._nbNodes;
2021 for (i = nbInitNod+1; i <= nbNod; ++i )
2023 const netgen::MeshPoint& ngPoint = ngMesh.Point(i);
2024 SMDS_MeshNode* node = NULL;
2025 TopoDS_Vertex aVert;
2026 // First, netgen creates nodes on vertices in occgeo.vmap,
2027 // so node index corresponds to vertex index
2028 // but (issue 0020776) netgen does not create nodes with equal coordinates
2029 if ( i-nbInitNod <= occgeo.vmap.Extent() )
2031 gp_Pnt p ( NGPOINT_COORDS(ngPoint) );
2032 for (int iV = i-nbInitNod; aVert.IsNull() && iV <= occgeo.vmap.Extent(); ++iV)
2034 aVert = TopoDS::Vertex( occgeo.vmap( iV ));
2035 gp_Pnt pV = BRep_Tool::Pnt( aVert );
2036 if ( p.SquareDistance( pV ) > 1e-20 )
2039 node = const_cast<SMDS_MeshNode*>( SMESH_Algo::VertexNode( aVert, meshDS ));
2042 if (!node) // node not found on vertex
2044 node = meshDS->AddNode( NGPOINT_COORDS( ngPoint ));
2045 if (!aVert.IsNull())
2046 meshDS->SetNodeOnVertex(node, aVert);
2051 // -------------------------------------------
2052 // Create mesh segments along geometric edges
2053 // -------------------------------------------
2055 int nbInitSeg = initState._nbSegments;
2056 for (i = nbInitSeg+1; i <= nbSeg; ++i )
2058 const netgen::Segment& seg = ngMesh.LineSegment(i);
2060 int pinds[3] = { seg.pnums[0], seg.pnums[1], seg.pnums[2] };
2063 for (int j=0; j < 3; ++j)
2065 int pind = pinds[j];
2066 if (pind <= 0 || !nodeVec_ACCESS(pind))
2074 int aGeomEdgeInd = seg.epgeominfo[j].edgenr;
2075 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
2076 aEdge = TopoDS::Edge(occgeo.emap(aGeomEdgeInd));
2078 param = seg.epgeominfo[j].dist;
2081 else // middle point
2083 param = param2 * 0.5;
2085 if (!aEdge.IsNull() && nodeVec_ACCESS(pind)->getshapeId() < 1)
2087 meshDS->SetNodeOnEdge(nodeVec_ACCESS(pind), aEdge, param);
2092 SMDS_MeshEdge* edge = 0;
2093 if (nbp == 2) // second order ?
2095 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1])))
2097 if ( quadHelper ) // final mesh must be quadratic
2098 edge = quadHelper->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
2100 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
2104 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
2105 nodeVec_ACCESS(pinds[2])))
2107 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
2108 nodeVec_ACCESS(pinds[2]));
2112 if ( comment.empty() ) comment << "Cannot create a mesh edge";
2113 MESSAGE("Cannot create a mesh edge");
2114 nbSeg = nbFac = nbVol = 0;
2117 if ( !aEdge.IsNull() && edge->getshapeId() < 1 )
2118 meshDS->SetMeshElementOnShape(edge, aEdge);
2120 else if ( comment.empty() )
2122 comment << "Invalid netgen segment #" << i;
2126 // ----------------------------------------
2127 // Create mesh faces along geometric faces
2128 // ----------------------------------------
2130 int nbInitFac = initState._nbFaces;
2131 int quadFaceID = ngMesh.GetNFD() + 1;
2132 if ( nbInitFac < nbFac )
2133 // add a faces descriptor to exclude qudrangle elements generated by NETGEN
2134 // from computation of 3D mesh
2135 ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(quadFaceID, /*solid1=*/0, /*solid2=*/0, 0));
2137 vector<const SMDS_MeshNode*> nodes;
2138 for (i = nbInitFac+1; i <= nbFac; ++i )
2140 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
2141 const int aGeomFaceInd = elem.GetIndex();
2143 if (aGeomFaceInd > 0 && aGeomFaceInd <= occgeo.fmap.Extent())
2144 aFace = TopoDS::Face(occgeo.fmap(aGeomFaceInd));
2146 for ( int j = 1; j <= elem.GetNP(); ++j )
2148 int pind = elem.PNum(j);
2149 if ( pind < 1 || pind >= (int) nodeVec.size() )
2151 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind))
2153 nodes.push_back( node );
2154 if (!aFace.IsNull() && node->getshapeId() < 1)
2156 const netgen::PointGeomInfo& pgi = elem.GeomInfoPi(j);
2157 meshDS->SetNodeOnFace(node, aFace, pgi.u, pgi.v);
2161 if ((int) nodes.size() != elem.GetNP() )
2163 if ( comment.empty() )
2164 comment << "Invalid netgen 2d element #" << i;
2165 continue; // bad node ids
2167 SMDS_MeshFace* face = NULL;
2168 switch (elem.GetType())
2171 if ( quadHelper ) // final mesh must be quadratic
2172 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2]);
2174 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2]);
2177 if ( quadHelper ) // final mesh must be quadratic
2178 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2180 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2181 // exclude qudrangle elements from computation of 3D mesh
2182 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2185 nodes[5] = mediumNode( nodes[0],nodes[1],nodes[5], quadHelper );
2186 nodes[3] = mediumNode( nodes[1],nodes[2],nodes[3], quadHelper );
2187 nodes[4] = mediumNode( nodes[2],nodes[0],nodes[4], quadHelper );
2188 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[5],nodes[3],nodes[4]);
2191 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2192 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2193 nodes[5] = mediumNode( nodes[2],nodes[3],nodes[5], quadHelper );
2194 nodes[6] = mediumNode( nodes[3],nodes[0],nodes[6], quadHelper );
2195 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3],
2196 nodes[4],nodes[7],nodes[5],nodes[6]);
2197 // exclude qudrangle elements from computation of 3D mesh
2198 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2201 MESSAGE("NETGEN created a face of unexpected type, ignoring");
2206 if ( comment.empty() ) comment << "Cannot create a mesh face";
2207 MESSAGE("Cannot create a mesh face");
2208 nbSeg = nbFac = nbVol = 0;
2211 if ( !aFace.IsNull() )
2212 meshDS->SetMeshElementOnShape( face, aFace );
2215 // ------------------
2216 // Create tetrahedra
2217 // ------------------
2219 for ( i = 1; i <= nbVol; ++i )
2221 const netgen::Element& elem = ngMesh.VolumeElement(i);
2222 int aSolidInd = elem.GetIndex();
2223 TopoDS_Solid aSolid;
2224 if ( aSolidInd > 0 && aSolidInd <= occgeo.somap.Extent() )
2225 aSolid = TopoDS::Solid(occgeo.somap(aSolidInd));
2227 for ( int j = 1; j <= elem.GetNP(); ++j )
2229 int pind = elem.PNum(j);
2230 if ( pind < 1 || pind >= (int)nodeVec.size() )
2232 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind) )
2234 nodes.push_back(node);
2235 if ( !aSolid.IsNull() && node->getshapeId() < 1 )
2236 meshDS->SetNodeInVolume(node, aSolid);
2239 if ((int) nodes.size() != elem.GetNP() )
2241 if ( comment.empty() )
2242 comment << "Invalid netgen 3d element #" << i;
2245 SMDS_MeshVolume* vol = NULL;
2246 switch ( elem.GetType() )
2249 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3]);
2252 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2253 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2254 nodes[5] = mediumNode( nodes[2],nodes[0],nodes[5], quadHelper );
2255 nodes[6] = mediumNode( nodes[0],nodes[3],nodes[6], quadHelper );
2256 nodes[8] = mediumNode( nodes[1],nodes[3],nodes[8], quadHelper );
2257 nodes[9] = mediumNode( nodes[2],nodes[3],nodes[9], quadHelper );
2258 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3],
2259 nodes[4],nodes[7],nodes[5],nodes[6],nodes[8],nodes[9]);
2262 MESSAGE("NETGEN created a volume of unexpected type, ignoring");
2267 if ( comment.empty() ) comment << "Cannot create a mesh volume";
2268 MESSAGE("Cannot create a mesh volume");
2269 nbSeg = nbFac = nbVol = 0;
2272 if (!aSolid.IsNull())
2273 meshDS->SetMeshElementOnShape(vol, aSolid);
2275 return comment.empty() ? 0 : 1;
2280 //================================================================================
2282 * \brief Restrict size of elements on the given edge
2284 //================================================================================
2286 void setLocalSize(const TopoDS_Edge& edge,
2290 if ( size <= std::numeric_limits<double>::min() )
2292 Standard_Real u1, u2;
2293 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, u1, u2);
2294 if ( curve.IsNull() )
2296 TopoDS_Iterator vIt( edge );
2297 if ( !vIt.More() ) return;
2298 gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex( vIt.Value() ));
2299 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
2303 const int nb = (int)( 1.5 * SMESH_Algo::EdgeLength( edge ) / size );
2304 Standard_Real delta = (u2-u1)/nb;
2305 for(int i=0; i<nb; i++)
2307 Standard_Real u = u1 + delta*i;
2308 gp_Pnt p = curve->Value(u);
2309 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
2310 netgen::Point3d pi(p.X(), p.Y(), p.Z());
2311 double resultSize = mesh.GetH(pi);
2312 if ( resultSize - size > 0.1*size )
2313 // netgen does restriction iff oldH/newH > 1.2 (localh.cpp:136)
2314 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), resultSize/1.201 );
2319 //================================================================================
2321 * \brief Convert error into text
2323 //================================================================================
2325 std::string text(int err)
2330 SMESH_Comment("Error in netgen::OCCGenerateMesh() at ") << netgen::multithread.task;
2333 //================================================================================
2335 * \brief Convert exception into text
2337 //================================================================================
2339 std::string text(Standard_Failure& ex)
2341 SMESH_Comment str("Exception in netgen::OCCGenerateMesh()");
2342 str << " at " << netgen::multithread.task
2343 << ": " << ex.DynamicType()->Name();
2344 if ( ex.GetMessageString() && strlen( ex.GetMessageString() ))
2345 str << ": " << ex.GetMessageString();
2348 //================================================================================
2350 * \brief Convert exception into text
2352 //================================================================================
2354 std::string text(netgen::NgException& ex)
2356 SMESH_Comment str("NgException");
2357 if ( strlen( netgen::multithread.task ) > 0 )
2358 str << " at " << netgen::multithread.task;
2359 str << ": " << ex.What();
2363 //================================================================================
2365 * \brief Looks for triangles lying on a SOLID
2367 //================================================================================
2369 bool hasBadElemOnSolid( const list<const SMDS_MeshElement*>& elems,
2370 SMESH_subMesh* solidSM )
2372 TopTools_IndexedMapOfShape solidSubs;
2373 TopExp::MapShapes( solidSM->GetSubShape(), solidSubs );
2374 SMESHDS_Mesh* mesh = solidSM->GetFather()->GetMeshDS();
2376 list<const SMDS_MeshElement*>::const_iterator e = elems.begin();
2377 for ( ; e != elems.end(); ++e )
2379 const SMDS_MeshElement* elem = *e;
2380 // if ( elem->GetType() != SMDSAbs_Face ) -- 23047
2382 int nbNodesOnSolid = 0, nbNodes = elem->NbNodes();
2383 SMDS_NodeIteratorPtr nIt = elem->nodeIterator();
2384 while ( nIt->more() )
2386 const SMDS_MeshNode* n = nIt->next();
2387 const TopoDS_Shape& s = mesh->IndexToShape( n->getshapeId() );
2388 nbNodesOnSolid += ( !s.IsNull() && solidSubs.Contains( s ));
2389 if ( nbNodesOnSolid > 2 ||
2390 nbNodesOnSolid == nbNodes)
2397 const double edgeMeshingTime = 0.001;
2398 const double faceMeshingTime = 0.019;
2399 const double edgeFaceMeshingTime = edgeMeshingTime + faceMeshingTime;
2400 const double faceOptimizTime = 0.06;
2401 const double voluMeshingTime = 0.15;
2402 const double volOptimizeTime = 0.77;
2405 //=============================================================================
2407 * Here we are going to use the NETGEN mesher
2409 //=============================================================================
2411 bool NETGENPlugin_Mesher::Compute()
2413 NETGENPlugin_NetgenLibWrapper ngLib;
2415 netgen::MeshingParameters& mparams = netgen::mparam;
2416 MESSAGE("Compute with:\n"
2417 " max size = " << mparams.maxh << "\n"
2418 " segments per edge = " << mparams.segmentsperedge);
2420 " growth rate = " << mparams.grading << "\n"
2421 " elements per radius = " << mparams.curvaturesafety << "\n"
2422 " second order = " << mparams.secondorder << "\n"
2423 " quad allowed = " << mparams.quad << "\n"
2424 " surface curvature = " << mparams.uselocalh << "\n"
2425 " fuse edges = " << netgen::merge_solids);
2427 SMESH_ComputeErrorPtr error = SMESH_ComputeError::New();
2428 SMESH_MesherHelper quadHelper( *_mesh );
2429 quadHelper.SetIsQuadratic( mparams.secondorder );
2431 static string debugFile = "/tmp/ngMesh.py"; /* to call toPython( _ngMesh, debugFile )
2432 while debugging netgen */
2433 // -------------------------
2434 // Prepare OCC geometry
2435 // -------------------------
2437 netgen::OCCGeometry occgeo;
2438 list< SMESH_subMesh* > meshedSM[3]; // for 0-2 dimensions
2439 NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
2440 PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
2443 _totalTime = edgeFaceMeshingTime;
2445 _totalTime += faceOptimizTime;
2447 _totalTime += voluMeshingTime + ( _optimize ? volOptimizeTime : 0 );
2448 double doneTime = 0;
2451 _curShapeIndex = -1;
2453 // -------------------------
2454 // Generate the mesh
2455 // -------------------------
2458 NETGENPlugin_ngMeshInfo initState; // it remembers size of ng mesh equal to size of Smesh
2460 SMESH_Comment comment;
2463 // vector of nodes in which node index == netgen ID
2464 vector< const SMDS_MeshNode* > nodeVec;
2472 // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
2473 mparams.uselocalh = false;
2474 mparams.grading = 0.8; // not limitited size growth
2476 if ( _simpleHyp->GetNumberOfSegments() )
2478 mparams.maxh = occgeo.boundingbox.Diam();
2481 mparams.maxh = _simpleHyp->GetLocalLength();
2484 if ( mparams.maxh == 0.0 )
2485 mparams.maxh = occgeo.boundingbox.Diam();
2486 if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
2487 mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
2489 // Local size on faces
2490 occgeo.face_maxh = mparams.maxh;
2492 // Let netgen create _ngMesh and calculate element size on not meshed shapes
2496 int startWith = netgen::MESHCONST_ANALYSE;
2497 int endWith = netgen::MESHCONST_ANALYSE;
2502 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2504 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2506 if(netgen::multithread.terminate)
2509 comment << text(err);
2511 catch (Standard_Failure& ex)
2513 comment << text(ex);
2515 err = 0; //- MESHCONST_ANALYSE isn't so important step
2518 ngLib.setMesh(( Ng_Mesh*) _ngMesh );
2520 _ngMesh->ClearFaceDescriptors(); // we make descriptors our-self
2524 // Pass 1D simple parameters to NETGEN
2525 // --------------------------------
2526 int nbSeg = _simpleHyp->GetNumberOfSegments();
2527 double segSize = _simpleHyp->GetLocalLength();
2528 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
2530 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
2532 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
2533 setLocalSize( e, segSize, *_ngMesh );
2536 else // if ( ! _simpleHyp )
2538 // Local size on vertices and edges
2539 // --------------------------------
2540 for(std::map<int,double>::const_iterator it=EdgeId2LocalSize.begin(); it!=EdgeId2LocalSize.end(); it++)
2542 int key = (*it).first;
2543 double hi = (*it).second;
2544 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2545 const TopoDS_Edge& e = TopoDS::Edge(shape);
2546 setLocalSize( e, hi, *_ngMesh );
2548 for(std::map<int,double>::const_iterator it=VertexId2LocalSize.begin(); it!=VertexId2LocalSize.end(); it++)
2550 int key = (*it).first;
2551 double hi = (*it).second;
2552 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2553 const TopoDS_Vertex& v = TopoDS::Vertex(shape);
2554 gp_Pnt p = BRep_Tool::Pnt(v);
2555 NETGENPlugin_Mesher::RestrictLocalSize( *_ngMesh, p.XYZ(), hi );
2557 for(map<int,double>::const_iterator it=FaceId2LocalSize.begin();
2558 it!=FaceId2LocalSize.end(); it++)
2560 int key = (*it).first;
2561 double val = (*it).second;
2562 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2563 int faceNgID = occgeo.fmap.FindIndex(shape);
2564 occgeo.SetFaceMaxH(faceNgID, val);
2565 for ( TopExp_Explorer edgeExp( shape, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
2566 setLocalSize( TopoDS::Edge( edgeExp.Current() ), val, *_ngMesh );
2570 // Precompute internal edges (issue 0020676) in order to
2571 // add mesh on them correctly (twice) to netgen mesh
2572 if ( !err && internals.hasInternalEdges() )
2574 // load internal shapes into OCCGeometry
2575 netgen::OCCGeometry intOccgeo;
2576 internals.getInternalEdges( intOccgeo.fmap, intOccgeo.emap, intOccgeo.vmap, meshedSM );
2577 intOccgeo.boundingbox = occgeo.boundingbox;
2578 intOccgeo.shape = occgeo.shape;
2579 intOccgeo.face_maxh.SetSize(intOccgeo.fmap.Extent());
2580 intOccgeo.face_maxh = netgen::mparam.maxh;
2581 netgen::Mesh *tmpNgMesh = NULL;
2585 // compute local H on internal shapes in the main mesh
2586 //OCCSetLocalMeshSize(intOccgeo, *_ngMesh); it deletes _ngMesh->localH
2588 // let netgen create a temporary mesh
2590 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2592 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2594 if(netgen::multithread.terminate)
2597 // copy LocalH from the main to temporary mesh
2598 initState.transferLocalH( _ngMesh, tmpNgMesh );
2600 // compute mesh on internal edges
2601 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2603 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2605 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2607 comment << text(err);
2609 catch (Standard_Failure& ex)
2611 comment << text(ex);
2614 initState.restoreLocalH( tmpNgMesh );
2616 // fill SMESH by netgen mesh
2617 vector< const SMDS_MeshNode* > tmpNodeVec;
2618 FillSMesh( intOccgeo, *tmpNgMesh, initState, *_mesh, tmpNodeVec, comment );
2619 err = ( err || !comment.empty() );
2621 nglib::Ng_DeleteMesh((nglib::Ng_Mesh*)tmpNgMesh);
2624 // Fill _ngMesh with nodes and segments of computed submeshes
2627 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_0D ]) &&
2628 FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_1D ], &quadHelper));
2630 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2635 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2640 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2642 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2644 if(netgen::multithread.terminate)
2647 comment << text(err);
2649 catch (Standard_Failure& ex)
2651 comment << text(ex);
2656 _ticTime = ( doneTime += edgeMeshingTime ) / _totalTime / _progressTic;
2658 mparams.uselocalh = true; // restore as it is used at surface optimization
2660 // ---------------------
2661 // compute surface mesh
2662 // ---------------------
2665 // Pass 2D simple parameters to NETGEN
2667 if ( double area = _simpleHyp->GetMaxElementArea() ) {
2669 mparams.maxh = sqrt(2. * area/sqrt(3.0));
2670 mparams.grading = 0.4; // moderate size growth
2673 // length from edges
2674 if ( _ngMesh->GetNSeg() ) {
2675 double edgeLength = 0;
2676 TopTools_MapOfShape visitedEdges;
2677 for ( TopExp_Explorer exp( _shape, TopAbs_EDGE ); exp.More(); exp.Next() )
2678 if( visitedEdges.Add(exp.Current()) )
2679 edgeLength += SMESH_Algo::EdgeLength( TopoDS::Edge( exp.Current() ));
2680 // we have to multiply length by 2 since for each TopoDS_Edge there
2681 // are double set of NETGEN edges, in other words, we have to
2682 // divide _ngMesh->GetNSeg() by 2.
2683 mparams.maxh = 2*edgeLength / _ngMesh->GetNSeg();
2686 mparams.maxh = 1000;
2688 mparams.grading = 0.2; // slow size growth
2690 mparams.quad = _simpleHyp->GetAllowQuadrangles();
2691 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2692 _ngMesh->SetGlobalH (mparams.maxh);
2693 netgen::Box<3> bb = occgeo.GetBoundingBox();
2694 bb.Increase (bb.Diam()/20);
2695 _ngMesh->SetLocalH (bb.PMin(), bb.PMax(), mparams.grading);
2698 // Care of vertices internal in faces (issue 0020676)
2699 if ( internals.hasInternalVertexInFace() )
2701 // store computed segments in SMESH in order not to create SMESH
2702 // edges for ng segments added by AddIntVerticesInFaces()
2703 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2704 // add segments to faces with internal vertices
2705 AddIntVerticesInFaces( occgeo, *_ngMesh, nodeVec, internals );
2706 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2709 // Build viscous layers
2710 if ( _isViscousLayers2D ||
2711 StdMeshers_ViscousLayers2D::HasProxyMesh( TopoDS::Face( occgeo.fmap(1) ), *_mesh ))
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 if ( viscousMesh->NbProxySubMeshes() == 0 )
2727 // exclude from computation ng segments built on EDGEs of F
2728 for (int i = 1; i <= _ngMesh->GetNSeg(); i++)
2730 netgen::Segment & seg = _ngMesh->LineSegment(i);
2731 if (seg.si == faceID)
2734 // add new segments to _ngMesh instead of excluded ones
2735 helper.SetSubShape( F );
2737 StdMeshers_FaceSide::GetFaceWires( F, *_mesh, /*skipMediumNodes=*/true,
2738 error, viscousMesh );
2739 error = AddSegmentsToMesh( *_ngMesh, occgeo, wires, helper, nodeVec );
2741 if ( !error ) error = SMESH_ComputeError::New();
2743 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2746 // Let netgen compute 2D mesh
2747 startWith = netgen::MESHCONST_MESHSURFACE;
2748 endWith = _optimize ? netgen::MESHCONST_OPTSURFACE : netgen::MESHCONST_MESHSURFACE;
2753 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2755 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2757 if(netgen::multithread.terminate)
2760 comment << text (err);
2762 catch (Standard_Failure& ex)
2764 comment << text(ex);
2765 //err = 1; -- try to make volumes anyway
2767 catch (netgen::NgException exc)
2769 comment << text(exc);
2770 //err = 1; -- try to make volumes anyway
2775 doneTime += faceMeshingTime + ( _optimize ? faceOptimizTime : 0 );
2776 _ticTime = doneTime / _totalTime / _progressTic;
2778 // ---------------------
2779 // generate volume mesh
2780 // ---------------------
2781 // Fill _ngMesh with nodes and faces of computed 2D submeshes
2782 if ( !err && _isVolume && ( !meshedSM[ MeshDim_2D ].empty() || mparams.quad ))
2784 // load SMESH with computed segments and faces
2785 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2787 // compute pyramids on quadrangles
2788 SMESH_ProxyMesh::Ptr proxyMesh;
2789 if ( _mesh->NbQuadrangles() > 0 )
2790 for ( int iS = 1; iS <= occgeo.somap.Extent(); ++iS )
2792 StdMeshers_QuadToTriaAdaptor* Adaptor = new StdMeshers_QuadToTriaAdaptor;
2793 proxyMesh.reset( Adaptor );
2795 int nbPyrams = _mesh->NbPyramids();
2796 Adaptor->Compute( *_mesh, occgeo.somap(iS) );
2797 if ( nbPyrams != _mesh->NbPyramids() )
2799 list< SMESH_subMesh* > quadFaceSM;
2800 for (TopExp_Explorer face(occgeo.somap(iS), TopAbs_FACE); face.More(); face.Next())
2801 if ( Adaptor->GetProxySubMesh( face.Current() ))
2803 quadFaceSM.push_back( _mesh->GetSubMesh( face.Current() ));
2804 meshedSM[ MeshDim_2D ].remove( quadFaceSM.back() );
2806 FillNgMesh(occgeo, *_ngMesh, nodeVec, quadFaceSM, &quadHelper, proxyMesh);
2809 // fill _ngMesh with faces of sub-meshes
2810 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_2D ], &quadHelper));
2811 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2812 //toPython( _ngMesh, "/tmp/ngPython.py");
2814 if (!err && _isVolume)
2816 // Pass 3D simple parameters to NETGEN
2817 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
2818 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
2820 if ( double vol = simple3d->GetMaxElementVolume() ) {
2822 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
2823 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2826 // length from faces
2827 mparams.maxh = _ngMesh->AverageH();
2829 _ngMesh->SetGlobalH (mparams.maxh);
2830 mparams.grading = 0.4;
2832 _ngMesh->CalcLocalH(mparams.grading);
2834 _ngMesh->CalcLocalH();
2837 // Care of vertices internal in solids and internal faces (issue 0020676)
2838 if ( internals.hasInternalVertexInSolid() || internals.hasInternalFaces() )
2840 // store computed faces in SMESH in order not to create SMESH
2841 // faces for ng faces added here
2842 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2843 // add ng faces to solids with internal vertices
2844 AddIntVerticesInSolids( occgeo, *_ngMesh, nodeVec, internals );
2845 // duplicate mesh faces on internal faces
2846 FixIntFaces( occgeo, *_ngMesh, internals );
2847 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2849 // Let netgen compute 3D mesh
2850 startWith = endWith = netgen::MESHCONST_MESHVOLUME;
2855 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2857 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2859 if(netgen::multithread.terminate)
2862 if ( comment.empty() ) // do not overwrite a previos error
2863 comment << text(err);
2865 catch (Standard_Failure& ex)
2867 if ( comment.empty() ) // do not overwrite a previos error
2868 comment << text(ex);
2871 catch (netgen::NgException exc)
2873 if ( comment.empty() ) // do not overwrite a previos error
2874 comment << text(exc);
2877 _ticTime = ( doneTime += voluMeshingTime ) / _totalTime / _progressTic;
2879 // Let netgen optimize 3D mesh
2880 if ( !err && _optimize )
2882 startWith = endWith = netgen::MESHCONST_OPTVOLUME;
2887 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2889 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2891 if(netgen::multithread.terminate)
2894 if ( comment.empty() ) // do not overwrite a previos error
2895 comment << text(err);
2897 catch (Standard_Failure& ex)
2899 if ( comment.empty() ) // do not overwrite a previos error
2900 comment << text(ex);
2902 catch (netgen::NgException exc)
2904 if ( comment.empty() ) // do not overwrite a previos error
2905 comment << text(exc);
2909 if (!err && mparams.secondorder > 0)
2914 if ( !meshedSM[ MeshDim_1D ].empty() )
2916 // remove segments not attached to geometry (IPAL0052479)
2917 for (int i = 1; i <= _ngMesh->GetNSeg(); ++i)
2919 const netgen::Segment & seg = _ngMesh->LineSegment (i);
2920 if ( seg.epgeominfo[ 0 ].edgenr == 0 )
2921 _ngMesh->DeleteSegment( i );
2923 _ngMesh->Compress();
2925 // convert to quadratic
2926 netgen::OCCRefinementSurfaces ref (occgeo);
2927 ref.MakeSecondOrder (*_ngMesh);
2929 // care of elements already loaded to SMESH
2930 // if ( initState._nbSegments > 0 )
2931 // makeQuadratic( occgeo.emap, _mesh );
2932 // if ( initState._nbFaces > 0 )
2933 // makeQuadratic( occgeo.fmap, _mesh );
2935 catch (Standard_Failure& ex)
2937 if ( comment.empty() ) // do not overwrite a previos error
2938 comment << "Exception in netgen at passing to 2nd order ";
2940 catch (netgen::NgException exc)
2942 if ( comment.empty() ) // do not overwrite a previos error
2943 comment << exc.What();
2948 _ticTime = 0.98 / _progressTic;
2950 int nbNod = _ngMesh->GetNP();
2951 int nbSeg = _ngMesh->GetNSeg();
2952 int nbFac = _ngMesh->GetNSE();
2953 int nbVol = _ngMesh->GetNE();
2954 bool isOK = ( !err && (_isVolume ? (nbVol > 0) : (nbFac > 0)) );
2956 // Feed back the SMESHDS with the generated Nodes and Elements
2957 if ( true /*isOK*/ ) // get whatever built
2959 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2961 if ( quadHelper.GetIsQuadratic() ) // remove free nodes
2962 for ( size_t i = 0; i < nodeVec.size(); ++i )
2963 if ( nodeVec[i] && nodeVec[i]->NbInverseElements() == 0 )
2964 _mesh->GetMeshDS()->RemoveFreeNode( nodeVec[i], 0, /*fromGroups=*/false );
2966 SMESH_ComputeErrorPtr readErr = ReadErrors(nodeVec);
2967 if ( readErr && !readErr->myBadElements.empty() )
2970 if ( !comment.empty() && !readErr->myComment.empty() ) comment += "\n";
2971 comment += readErr->myComment;
2973 if ( error->IsOK() && ( !isOK || comment.size() > 0 ))
2974 error->myName = COMPERR_ALGO_FAILED;
2975 if ( !comment.empty() )
2976 error->myComment = comment;
2978 // SetIsAlwaysComputed( true ) to empty sub-meshes, which
2979 // appear if the geometry contains coincident sub-shape due
2980 // to bool merge_solids = 1; in netgen/libsrc/occ/occgenmesh.cpp
2981 const int nbMaps = 2;
2982 const TopTools_IndexedMapOfShape* geoMaps[nbMaps] =
2983 { & occgeo.vmap, & occgeo.emap/*, & occgeo.fmap*/ };
2984 for ( int iMap = 0; iMap < nbMaps; ++iMap )
2985 for (int i = 1; i <= geoMaps[iMap]->Extent(); i++)
2986 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( geoMaps[iMap]->FindKey(i)))
2987 if ( !sm->IsMeshComputed() )
2988 sm->SetIsAlwaysComputed( true );
2990 // set bad compute error to subshapes of all failed sub-shapes
2991 if ( !error->IsOK() )
2993 bool pb2D = false, pb3D = false;
2994 for (int i = 1; i <= occgeo.fmap.Extent(); i++) {
2995 int status = occgeo.facemeshstatus[i-1];
2996 if (status == netgen::FACE_MESHED_OK ) continue;
2997 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.fmap( i ))) {
2998 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
2999 if ( !smError || smError->IsOK() ) {
3000 if ( status == netgen::FACE_FAILED )
3001 smError.reset( new SMESH_ComputeError( *error ));
3003 smError.reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED, "Ignored" ));
3004 if ( SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
3005 smError->myName = COMPERR_WARNING;
3007 pb2D = pb2D || smError->IsKO();
3010 if ( !pb2D ) // all faces are OK
3011 for (int i = 1; i <= occgeo.somap.Extent(); i++)
3012 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.somap( i )))
3014 bool smComputed = nbVol && !sm->IsEmpty();
3015 if ( smComputed && internals.hasInternalVertexInSolid( sm->GetId() ))
3017 int nbIntV = internals.getSolidsWithVertices().find( sm->GetId() )->second.size();
3018 SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
3019 smComputed = ( smDS->NbElements() > 0 || smDS->NbNodes() > nbIntV );
3021 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
3022 if ( !smComputed && ( !smError || smError->IsOK() ))
3024 smError.reset( new SMESH_ComputeError( *error ));
3025 if ( nbVol && SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
3027 smError->myName = COMPERR_WARNING;
3029 else if ( !smError->myBadElements.empty() ) // bad surface mesh
3031 if ( !hasBadElemOnSolid( smError->myBadElements, sm ))
3035 pb3D = pb3D || ( smError && smError->IsKO() );
3037 if ( !pb2D && !pb3D )
3038 err = 0; // no fatal errors, only warnings
3041 ngLib._isComputeOk = !err;
3046 //=============================================================================
3050 //=============================================================================
3051 bool NETGENPlugin_Mesher::Evaluate(MapShapeNbElems& aResMap)
3053 netgen::MeshingParameters& mparams = netgen::mparam;
3056 // -------------------------
3057 // Prepare OCC geometry
3058 // -------------------------
3059 netgen::OCCGeometry occgeo;
3060 list< SMESH_subMesh* > meshedSM[4]; // for 0-3 dimensions
3061 NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
3062 PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
3064 bool tooManyElems = false;
3065 const int hugeNb = std::numeric_limits<int>::max() / 100;
3070 // pass 1D simple parameters to NETGEN
3073 // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
3074 mparams.uselocalh = false;
3075 mparams.grading = 0.8; // not limitited size growth
3077 if ( _simpleHyp->GetNumberOfSegments() )
3079 mparams.maxh = occgeo.boundingbox.Diam();
3082 mparams.maxh = _simpleHyp->GetLocalLength();
3085 if ( mparams.maxh == 0.0 )
3086 mparams.maxh = occgeo.boundingbox.Diam();
3087 if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
3088 mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
3090 // let netgen create _ngMesh and calculate element size on not meshed shapes
3091 NETGENPlugin_NetgenLibWrapper ngLib;
3092 netgen::Mesh *ngMesh = NULL;
3096 int startWith = netgen::MESHCONST_ANALYSE;
3097 int endWith = netgen::MESHCONST_MESHEDGES;
3099 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, mparams, startWith, endWith);
3101 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, startWith, endWith, optstr);
3104 if(netgen::multithread.terminate)
3107 ngLib.setMesh(( Ng_Mesh*) ngMesh );
3109 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( _shape ))
3110 sm->GetComputeError().reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED ));
3115 // Pass 1D simple parameters to NETGEN
3116 // --------------------------------
3117 int nbSeg = _simpleHyp->GetNumberOfSegments();
3118 double segSize = _simpleHyp->GetLocalLength();
3119 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
3121 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
3123 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
3124 setLocalSize( e, segSize, *ngMesh );
3127 else // if ( ! _simpleHyp )
3129 // Local size on vertices and edges
3130 // --------------------------------
3131 for(std::map<int,double>::const_iterator it=EdgeId2LocalSize.begin(); it!=EdgeId2LocalSize.end(); it++)
3133 int key = (*it).first;
3134 double hi = (*it).second;
3135 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3136 const TopoDS_Edge& e = TopoDS::Edge(shape);
3137 setLocalSize( e, hi, *ngMesh );
3139 for(std::map<int,double>::const_iterator it=VertexId2LocalSize.begin(); it!=VertexId2LocalSize.end(); it++)
3141 int key = (*it).first;
3142 double hi = (*it).second;
3143 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3144 const TopoDS_Vertex& v = TopoDS::Vertex(shape);
3145 gp_Pnt p = BRep_Tool::Pnt(v);
3146 NETGENPlugin_Mesher::RestrictLocalSize( *ngMesh, p.XYZ(), hi );
3148 for(map<int,double>::const_iterator it=FaceId2LocalSize.begin();
3149 it!=FaceId2LocalSize.end(); it++)
3151 int key = (*it).first;
3152 double val = (*it).second;
3153 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3154 int faceNgID = occgeo.fmap.FindIndex(shape);
3155 occgeo.SetFaceMaxH(faceNgID, val);
3156 for ( TopExp_Explorer edgeExp( shape, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
3157 setLocalSize( TopoDS::Edge( edgeExp.Current() ), val, *ngMesh );
3160 // calculate total nb of segments and length of edges
3161 double fullLen = 0.0;
3163 int entity = mparams.secondorder > 0 ? SMDSEntity_Quad_Edge : SMDSEntity_Edge;
3164 TopTools_DataMapOfShapeInteger Edge2NbSeg;
3165 for (TopExp_Explorer exp(_shape, TopAbs_EDGE); exp.More(); exp.Next())
3167 TopoDS_Edge E = TopoDS::Edge( exp.Current() );
3168 if( !Edge2NbSeg.Bind(E,0) )
3171 double aLen = SMESH_Algo::EdgeLength(E);
3174 vector<int>& aVec = aResMap[_mesh->GetSubMesh(E)];
3176 aVec.resize( SMDSEntity_Last, 0);
3178 fullNbSeg += aVec[ entity ];
3181 // store nb of segments computed by Netgen
3182 NCollection_Map<Link> linkMap;
3183 for (int i = 1; i <= ngMesh->GetNSeg(); ++i )
3185 const netgen::Segment& seg = ngMesh->LineSegment(i);
3186 Link link(seg[0], seg[1]);
3187 if ( !linkMap.Add( link )) continue;
3188 int aGeomEdgeInd = seg.epgeominfo[0].edgenr;
3189 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
3191 vector<int>& aVec = aResMap[_mesh->GetSubMesh(occgeo.emap(aGeomEdgeInd))];
3195 // store nb of nodes on edges computed by Netgen
3196 TopTools_DataMapIteratorOfDataMapOfShapeInteger Edge2NbSegIt(Edge2NbSeg);
3197 for (; Edge2NbSegIt.More(); Edge2NbSegIt.Next())
3199 vector<int>& aVec = aResMap[_mesh->GetSubMesh(Edge2NbSegIt.Key())];
3200 if ( aVec[ entity ] > 1 && aVec[ SMDSEntity_Node ] == 0 )
3201 aVec[SMDSEntity_Node] = mparams.secondorder > 0 ? 2*aVec[ entity ]-1 : aVec[ entity ]-1;
3203 fullNbSeg += aVec[ entity ];
3204 Edge2NbSeg( Edge2NbSegIt.Key() ) = aVec[ entity ];
3206 if ( fullNbSeg == 0 )
3213 if ( double area = _simpleHyp->GetMaxElementArea() ) {
3215 mparams.maxh = sqrt(2. * area/sqrt(3.0));
3216 mparams.grading = 0.4; // moderate size growth
3219 // length from edges
3220 mparams.maxh = fullLen/fullNbSeg;
3221 mparams.grading = 0.2; // slow size growth
3224 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3225 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3227 for (TopExp_Explorer exp(_shape, TopAbs_FACE); exp.More(); exp.Next())
3229 TopoDS_Face F = TopoDS::Face( exp.Current() );
3230 SMESH_subMesh *sm = _mesh->GetSubMesh(F);
3232 BRepGProp::SurfaceProperties(F,G);
3233 double anArea = G.Mass();
3234 tooManyElems = tooManyElems || ( anArea/hugeNb > mparams.maxh*mparams.maxh );
3236 if ( !tooManyElems )
3238 TopTools_MapOfShape egdes;
3239 for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next())
3240 if ( egdes.Add( exp1.Current() ))
3241 nb1d += Edge2NbSeg.Find(exp1.Current());
3243 int nbFaces = tooManyElems ? hugeNb : int( 4*anArea / (mparams.maxh*mparams.maxh*sqrt(3.)));
3244 int nbNodes = tooManyElems ? hugeNb : (( nbFaces*3 - (nb1d-1)*2 ) / 6 + 1 );
3246 vector<int> aVec(SMDSEntity_Last, 0);
3247 if( mparams.secondorder > 0 ) {
3248 int nb1d_in = (nbFaces*3 - nb1d) / 2;
3249 aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
3250 aVec[SMDSEntity_Quad_Triangle] = nbFaces;
3253 aVec[SMDSEntity_Node] = Max ( nbNodes, 0 );
3254 aVec[SMDSEntity_Triangle] = nbFaces;
3256 aResMap[sm].swap(aVec);
3263 // pass 3D simple parameters to NETGEN
3264 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
3265 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
3267 if ( double vol = simple3d->GetMaxElementVolume() ) {
3269 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
3270 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3273 // using previous length from faces
3275 mparams.grading = 0.4;
3276 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3279 BRepGProp::VolumeProperties(_shape,G);
3280 double aVolume = G.Mass();
3281 double tetrVol = 0.1179*mparams.maxh*mparams.maxh*mparams.maxh;
3282 tooManyElems = tooManyElems || ( aVolume/hugeNb > tetrVol );
3283 int nbVols = tooManyElems ? hugeNb : int(aVolume/tetrVol);
3284 int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
3285 vector<int> aVec(SMDSEntity_Last, 0 );
3286 if ( tooManyElems ) // avoid FPE
3288 aVec[SMDSEntity_Node] = hugeNb;
3289 aVec[ mparams.secondorder > 0 ? SMDSEntity_Quad_Tetra : SMDSEntity_Tetra] = hugeNb;
3293 if( mparams.secondorder > 0 ) {
3294 aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
3295 aVec[SMDSEntity_Quad_Tetra] = nbVols;
3298 aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
3299 aVec[SMDSEntity_Tetra] = nbVols;
3302 SMESH_subMesh *sm = _mesh->GetSubMesh(_shape);
3303 aResMap[sm].swap(aVec);
3309 double NETGENPlugin_Mesher::GetProgress(const SMESH_Algo* holder,
3310 const int * algoProgressTic,
3311 const double * algoProgress) const
3313 ((int&) _progressTic ) = *algoProgressTic + 1;
3315 if ( !_occgeom ) return 0;
3317 double progress = -1;
3320 if ( _ticTime < 0 && netgen::multithread.task[0] == 'O'/*Optimizing surface*/ )
3322 ((double&) _ticTime ) = edgeFaceMeshingTime / _totalTime / _progressTic;
3324 else if ( !_optimize /*&& _occgeom->fmap.Extent() > 1*/ )
3326 int doneShapeIndex = -1;
3327 while ( doneShapeIndex+1 < _occgeom->facemeshstatus.Size() &&
3328 _occgeom->facemeshstatus[ doneShapeIndex+1 ])
3330 if ( doneShapeIndex+1 != _curShapeIndex )
3332 ((int&) _curShapeIndex) = doneShapeIndex+1;
3333 double doneShapeRate = _curShapeIndex / double( _occgeom->fmap.Extent() );
3334 double doneTime = edgeMeshingTime + doneShapeRate * faceMeshingTime;
3335 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3336 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3337 // << " " << doneTime / _totalTime / _progressTic << endl;
3341 else if ( !_optimize && _occgeom->somap.Extent() > 1 )
3343 int curShapeIndex = _curShapeIndex;
3344 if ( _ngMesh->GetNE() > 0 )
3346 netgen::Element el = (*_ngMesh)[netgen::ElementIndex( _ngMesh->GetNE()-1 )];
3347 curShapeIndex = el.GetIndex();
3349 if ( curShapeIndex != _curShapeIndex )
3351 ((int&) _curShapeIndex) = curShapeIndex;
3352 double doneShapeRate = _curShapeIndex / double( _occgeom->somap.Extent() );
3353 double doneTime = edgeFaceMeshingTime + doneShapeRate * voluMeshingTime;
3354 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3355 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3356 // << " " << doneTime / _totalTime / _progressTic << endl;
3360 progress = Max( *algoProgressTic * _ticTime, *algoProgress );
3363 ((int&) *algoProgressTic )++;
3364 ((double&) *algoProgress) = progress;
3366 //cout << progress << " " << *algoProgressTic << " " << netgen::multithread.task << " "<< _ticTime << endl;
3368 return Min( progress, 0.99 );
3371 //================================================================================
3373 * \brief Read mesh entities preventing successful computation from "test.out" file
3375 //================================================================================
3377 SMESH_ComputeErrorPtr
3378 NETGENPlugin_Mesher::ReadErrors(const vector<const SMDS_MeshNode* >& nodeVec)
3380 SMESH_ComputeErrorPtr err = SMESH_ComputeError::New
3381 (COMPERR_BAD_INPUT_MESH, "Some edges multiple times in surface mesh");
3382 SMESH_File file("test.out");
3384 vector<int> three1(3), three2(3);
3385 const char* badEdgeStr = " multiple times in surface mesh";
3386 const int badEdgeStrLen = strlen( badEdgeStr );
3387 const int nbNodes = nodeVec.size();
3389 while( !file.eof() )
3391 if ( strncmp( file, "Edge ", 5 ) == 0 &&
3392 file.getInts( two ) &&
3393 strncmp( file, badEdgeStr, badEdgeStrLen ) == 0 &&
3394 two[0] < nbNodes && two[1] < nbNodes )
3396 err->myBadElements.push_back( new SMDS_LinearEdge( nodeVec[ two[0]], nodeVec[ two[1]] ));
3397 file += badEdgeStrLen;
3399 else if ( strncmp( file, "Intersecting: ", 14 ) == 0 )
3402 // openelement 18 with open element 126
3406 const char* pos = file;
3407 bool ok = ( strncmp( file, "openelement ", 12 ) == 0 );
3408 ok = ok && file.getInts( two );
3409 ok = ok && file.getInts( three1 );
3410 ok = ok && file.getInts( three2 );
3411 for ( int i = 0; ok && i < 3; ++i )
3412 ok = ( three1[i] < nbNodes && nodeVec[ three1[i]]);
3413 for ( int i = 0; ok && i < 3; ++i )
3414 ok = ( three2[i] < nbNodes && nodeVec[ three2[i]]);
3417 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three1[0]],
3418 nodeVec[ three1[1]],
3419 nodeVec[ three1[2]]));
3420 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three2[0]],
3421 nodeVec[ three2[1]],
3422 nodeVec[ three2[2]]));
3423 err->myComment = "Intersecting triangles";
3437 size_t nbBadElems = err->myBadElements.size();
3444 //================================================================================
3446 * \brief Write a python script creating an equivalent SALOME mesh.
3447 * This is useful to see what mesh is passed as input for the next step of mesh
3448 * generation (of mesh of higher dimension)
3450 //================================================================================
3452 void NETGENPlugin_Mesher::toPython( const netgen::Mesh* ngMesh,
3453 const std::string& pyFile)
3455 ofstream outfile(pyFile.c_str(), ios::out);
3456 if ( !outfile ) return;
3458 outfile << "import SMESH" << endl
3459 << "from salome.smesh import smeshBuilder" << endl
3460 << "smesh = smeshBuilder.New(salome.myStudy)" << endl
3461 << "mesh = smesh.Mesh()" << endl << endl;
3463 using namespace netgen;
3465 for (pi = PointIndex::BASE;
3466 pi < ngMesh->GetNP()+PointIndex::BASE; pi++)
3468 outfile << "mesh.AddNode( ";
3469 outfile << (*ngMesh)[pi](0) << ", ";
3470 outfile << (*ngMesh)[pi](1) << ", ";
3471 outfile << (*ngMesh)[pi](2) << ") ## "<< pi << endl;
3474 int nbDom = ngMesh->GetNDomains();
3475 for ( int i = 0; i < nbDom; ++i )
3476 outfile<< "grp" << i+1 << " = mesh.CreateEmptyGroup( SMESH.FACE, 'domain"<< i+1 << "')"<< endl;
3478 SurfaceElementIndex sei;
3479 for (sei = 0; sei < ngMesh->GetNSE(); sei++)
3481 outfile << "mesh.AddFace([ ";
3482 Element2d sel = (*ngMesh)[sei];
3483 for (int j = 0; j < sel.GetNP(); j++)
3484 outfile << sel[j] << ( j+1 < sel.GetNP() ? ", " : " ])");
3485 if ( sel.IsDeleted() ) outfile << " ## IsDeleted ";
3488 if ((*ngMesh)[sei].GetIndex())
3490 if ( int dom1 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainIn())
3491 outfile << "grp"<< dom1 <<".Add([ " << (int)sei+1 << " ])" << endl;
3492 if ( int dom2 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainOut())
3493 outfile << "grp"<< dom2 <<".Add([ " << (int)sei+1 << " ])" << endl;
3497 for (ElementIndex ei = 0; ei < ngMesh->GetNE(); ei++)
3499 Element el = (*ngMesh)[ei];
3500 outfile << "mesh.AddVolume([ ";
3501 for (int j = 0; j < el.GetNP(); j++)
3502 outfile << el[j] << ( j+1 < el.GetNP() ? ", " : " ])");
3506 for (int i = 1; i <= ngMesh->GetNSeg(); i++)
3508 const Segment & seg = ngMesh->LineSegment (i);
3509 outfile << "mesh.AddEdge([ "
3511 << seg[1] << " ])" << endl;
3513 cout << "Write " << pyFile << endl;
3516 //================================================================================
3518 * \brief Constructor of NETGENPlugin_ngMeshInfo
3520 //================================================================================
3522 NETGENPlugin_ngMeshInfo::NETGENPlugin_ngMeshInfo( netgen::Mesh* ngMesh):
3527 _nbNodes = ngMesh->GetNP();
3528 _nbSegments = ngMesh->GetNSeg();
3529 _nbFaces = ngMesh->GetNSE();
3530 _nbVolumes = ngMesh->GetNE();
3534 _nbNodes = _nbSegments = _nbFaces = _nbVolumes = 0;
3538 //================================================================================
3540 * \brief Copy LocalH member from one netgen mesh to another
3542 //================================================================================
3544 void NETGENPlugin_ngMeshInfo::transferLocalH( netgen::Mesh* fromMesh,
3545 netgen::Mesh* toMesh )
3547 if ( !fromMesh->LocalHFunctionGenerated() ) return;
3548 if ( !toMesh->LocalHFunctionGenerated() )
3550 toMesh->CalcLocalH(netgen::mparam.grading);
3552 toMesh->CalcLocalH();
3555 const size_t size = sizeof( netgen::LocalH );
3556 _copyOfLocalH = new char[ size ];
3557 memcpy( (void*)_copyOfLocalH, (void*)&toMesh->LocalHFunction(), size );
3558 memcpy( (void*)&toMesh->LocalHFunction(), (void*)&fromMesh->LocalHFunction(), size );
3561 //================================================================================
3563 * \brief Restore LocalH member of a netgen mesh
3565 //================================================================================
3567 void NETGENPlugin_ngMeshInfo::restoreLocalH( netgen::Mesh* toMesh )
3569 if ( _copyOfLocalH )
3571 const size_t size = sizeof( netgen::LocalH );
3572 memcpy( (void*)&toMesh->LocalHFunction(), (void*)_copyOfLocalH, size );
3573 delete [] _copyOfLocalH;
3578 //================================================================================
3580 * \brief Find "internal" sub-shapes
3582 //================================================================================
3584 NETGENPlugin_Internals::NETGENPlugin_Internals( SMESH_Mesh& mesh,
3585 const TopoDS_Shape& shape,
3587 : _mesh( mesh ), _is3D( is3D )
3589 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3591 TopExp_Explorer f,e;
3592 for ( f.Init( shape, TopAbs_FACE ); f.More(); f.Next() )
3594 int faceID = meshDS->ShapeToIndex( f.Current() );
3596 // find not computed internal edges
3598 for ( e.Init( f.Current().Oriented(TopAbs_FORWARD), TopAbs_EDGE ); e.More(); e.Next() )
3599 if ( e.Current().Orientation() == TopAbs_INTERNAL )
3601 SMESH_subMesh* eSM = mesh.GetSubMesh( e.Current() );
3602 if ( eSM->IsEmpty() )
3604 _e2face.insert( make_pair( eSM->GetId(), faceID ));
3605 for ( TopoDS_Iterator v(e.Current()); v.More(); v.Next() )
3606 _e2face.insert( make_pair( meshDS->ShapeToIndex( v.Value() ), faceID ));
3610 // find internal vertices in a face
3611 set<int> intVV; // issue 0020850 where same vertex is twice in a face
3612 for ( TopoDS_Iterator fSub( f.Current() ); fSub.More(); fSub.Next())
3613 if ( fSub.Value().ShapeType() == TopAbs_VERTEX )
3615 int vID = meshDS->ShapeToIndex( fSub.Value() );
3616 if ( intVV.insert( vID ).second )
3617 _f2v[ faceID ].push_back( vID );
3622 // find internal faces and their subshapes where nodes are to be doubled
3623 // to make a crack with non-sewed borders
3625 if ( f.Current().Orientation() == TopAbs_INTERNAL )
3627 _intShapes.insert( meshDS->ShapeToIndex( f.Current() ));
3630 list< TopoDS_Shape > edges;
3631 for ( e.Init( f.Current(), TopAbs_EDGE ); e.More(); e.Next())
3632 if ( SMESH_MesherHelper::NbAncestors( e.Current(), mesh, TopAbs_FACE ) > 1 )
3634 _intShapes.insert( meshDS->ShapeToIndex( e.Current() ));
3635 edges.push_back( e.Current() );
3636 // find border faces
3637 PShapeIteratorPtr fIt =
3638 SMESH_MesherHelper::GetAncestors( edges.back(),mesh,TopAbs_FACE );
3639 while ( const TopoDS_Shape* pFace = fIt->next() )
3640 if ( !pFace->IsSame( f.Current() ))
3641 _borderFaces.insert( meshDS->ShapeToIndex( *pFace ));
3644 // we consider vertex internal if it is shared by more than one internal edge
3645 list< TopoDS_Shape >::iterator edge = edges.begin();
3646 for ( ; edge != edges.end(); ++edge )
3647 for ( TopoDS_Iterator v( *edge ); v.More(); v.Next() )
3649 set<int> internalEdges;
3650 PShapeIteratorPtr eIt =
3651 SMESH_MesherHelper::GetAncestors( v.Value(),mesh,TopAbs_EDGE );
3652 while ( const TopoDS_Shape* pEdge = eIt->next() )
3654 int edgeID = meshDS->ShapeToIndex( *pEdge );
3655 if ( isInternalShape( edgeID ))
3656 internalEdges.insert( edgeID );
3658 if ( internalEdges.size() > 1 )
3659 _intShapes.insert( meshDS->ShapeToIndex( v.Value() ));
3663 } // loop on geom faces
3665 // find vertices internal in solids
3668 for ( TopExp_Explorer so(shape, TopAbs_SOLID); so.More(); so.Next())
3670 int soID = meshDS->ShapeToIndex( so.Current() );
3671 for ( TopoDS_Iterator soSub( so.Current() ); soSub.More(); soSub.Next())
3672 if ( soSub.Value().ShapeType() == TopAbs_VERTEX )
3673 _s2v[ soID ].push_back( meshDS->ShapeToIndex( soSub.Value() ));
3678 //================================================================================
3680 * \brief Find mesh faces on non-internal geom faces sharing internal edge
3681 * some nodes of which are to be doubled to make the second border of the "crack"
3683 //================================================================================
3685 void NETGENPlugin_Internals::findBorderElements( TIDSortedElemSet & borderElems )
3687 if ( _intShapes.empty() ) return;
3689 SMESH_Mesh& mesh = const_cast<SMESH_Mesh&>(_mesh);
3690 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3692 // loop on internal geom edges
3693 set<int>::const_iterator intShapeId = _intShapes.begin();
3694 for ( ; intShapeId != _intShapes.end(); ++intShapeId )
3696 const TopoDS_Shape& s = meshDS->IndexToShape( *intShapeId );
3697 if ( s.ShapeType() != TopAbs_EDGE ) continue;
3699 // get internal and non-internal geom faces sharing the internal edge <s>
3701 set<int>::iterator bordFace = _borderFaces.end();
3702 PShapeIteratorPtr faces = SMESH_MesherHelper::GetAncestors( s, _mesh, TopAbs_FACE );
3703 while ( const TopoDS_Shape* pFace = faces->next() )
3705 int faceID = meshDS->ShapeToIndex( *pFace );
3706 if ( isInternalShape( faceID ))
3709 bordFace = _borderFaces.insert( faceID ).first;
3711 if ( bordFace == _borderFaces.end() || !intFace ) continue;
3713 // get all links of mesh faces on internal geom face sharing nodes on edge <s>
3714 set< SMESH_OrientedLink > links; //!< links of faces on internal geom face
3715 list<const SMDS_MeshElement*> suspectFaces[2]; //!< mesh faces on border geom faces
3716 int nbSuspectFaces = 0;
3717 SMESHDS_SubMesh* intFaceSM = meshDS->MeshElements( intFace );
3718 if ( !intFaceSM || intFaceSM->NbElements() == 0 ) continue;
3719 SMESH_subMeshIteratorPtr smIt = mesh.GetSubMesh( s )->getDependsOnIterator(true,true);
3720 while ( smIt->more() )
3722 SMESHDS_SubMesh* sm = smIt->next()->GetSubMeshDS();
3723 if ( !sm ) continue;
3724 SMDS_NodeIteratorPtr nIt = sm->GetNodes();
3725 while ( nIt->more() )
3727 const SMDS_MeshNode* nOnEdge = nIt->next();
3728 SMDS_ElemIteratorPtr fIt = nOnEdge->GetInverseElementIterator(SMDSAbs_Face);
3729 while ( fIt->more() )
3731 const SMDS_MeshElement* f = fIt->next();
3732 const int nbNodes = f->NbCornerNodes();
3733 if ( intFaceSM->Contains( f ))
3735 for ( int i = 0; i < nbNodes; ++i )
3736 links.insert( SMESH_OrientedLink( f->GetNode(i), f->GetNode((i+1)%nbNodes)));
3741 for ( int i = 0; i < nbNodes; ++i )
3742 nbDblNodes += isInternalShape( f->GetNode(i)->getshapeId() );
3744 suspectFaces[ nbDblNodes < 2 ].push_back( f );
3750 // suspectFaces[0] having link with same orientation as mesh faces on
3751 // the internal geom face are <borderElems>. suspectFaces[1] have
3752 // only one node on edge <s>, we decide on them later (at the 2nd loop)
3753 // by links of <borderElems> found at the 1st and 2nd loops
3754 set< SMESH_OrientedLink > borderLinks;
3755 for ( int isPostponed = 0; isPostponed < 2; ++isPostponed )
3757 list<const SMDS_MeshElement*>::iterator fIt = suspectFaces[isPostponed].begin();
3758 for ( int nbF = 0; fIt != suspectFaces[isPostponed].end(); ++fIt, ++nbF )
3760 const SMDS_MeshElement* f = *fIt;
3761 bool isBorder = false, linkFound = false, borderLinkFound = false;
3762 list< SMESH_OrientedLink > faceLinks;
3763 int nbNodes = f->NbCornerNodes();
3764 for ( int i = 0; i < nbNodes; ++i )
3766 SMESH_OrientedLink link( f->GetNode(i), f->GetNode((i+1)%nbNodes));
3767 faceLinks.push_back( link );
3770 set< SMESH_OrientedLink >::iterator foundLink = links.find( link );
3771 if ( foundLink != links.end() )
3774 isBorder = ( foundLink->_reversed == link._reversed );
3775 if ( !isBorder && !isPostponed ) break;
3776 faceLinks.pop_back();
3778 else if ( isPostponed && !borderLinkFound )
3780 foundLink = borderLinks.find( link );
3781 if ( foundLink != borderLinks.end() )
3783 borderLinkFound = true;
3784 isBorder = ( foundLink->_reversed != link._reversed );
3791 borderElems.insert( f );
3792 borderLinks.insert( faceLinks.begin(), faceLinks.end() );
3794 else if ( !linkFound && !borderLinkFound )
3796 suspectFaces[1].push_back( f );
3797 if ( nbF > 2 * nbSuspectFaces )
3798 break; // dead loop protection
3805 //================================================================================
3807 * \brief put internal shapes in maps and fill in submeshes to precompute
3809 //================================================================================
3811 void NETGENPlugin_Internals::getInternalEdges( TopTools_IndexedMapOfShape& fmap,
3812 TopTools_IndexedMapOfShape& emap,
3813 TopTools_IndexedMapOfShape& vmap,
3814 list< SMESH_subMesh* > smToPrecompute[])
3816 if ( !hasInternalEdges() ) return;
3817 map<int,int>::const_iterator ev_face = _e2face.begin();
3818 for ( ; ev_face != _e2face.end(); ++ev_face )
3820 const TopoDS_Shape& ev = _mesh.GetMeshDS()->IndexToShape( ev_face->first );
3821 const TopoDS_Shape& face = _mesh.GetMeshDS()->IndexToShape( ev_face->second );
3823 ( ev.ShapeType() == TopAbs_EDGE ? emap : vmap ).Add( ev );
3825 //cout<<"INTERNAL EDGE or VERTEX "<<ev_face->first<<" on face "<<ev_face->second<<endl;
3827 smToPrecompute[ MeshDim_1D ].push_back( _mesh.GetSubMeshContaining( ev_face->first ));
3831 //================================================================================
3833 * \brief return shapes and submeshes to be meshed and already meshed boundary submeshes
3835 //================================================================================
3837 void NETGENPlugin_Internals::getInternalFaces( TopTools_IndexedMapOfShape& fmap,
3838 TopTools_IndexedMapOfShape& emap,
3839 list< SMESH_subMesh* >& intFaceSM,
3840 list< SMESH_subMesh* >& boundarySM)
3842 if ( !hasInternalFaces() ) return;
3844 // <fmap> and <emap> are for not yet meshed shapes
3845 // <intFaceSM> is for submeshes of faces
3846 // <boundarySM> is for meshed edges and vertices
3851 set<int> shapeIDs ( _intShapes );
3852 if ( !_borderFaces.empty() )
3853 shapeIDs.insert( _borderFaces.begin(), _borderFaces.end() );
3855 set<int>::const_iterator intS = shapeIDs.begin();
3856 for ( ; intS != shapeIDs.end(); ++intS )
3858 SMESH_subMesh* sm = _mesh.GetSubMeshContaining( *intS );
3860 if ( sm->GetSubShape().ShapeType() != TopAbs_FACE ) continue;
3862 intFaceSM.push_back( sm );
3864 // add submeshes of not computed internal faces
3865 if ( !sm->IsEmpty() ) continue;
3867 SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(true,true);
3868 while ( smIt->more() )
3871 const TopoDS_Shape& s = sm->GetSubShape();
3873 if ( sm->IsEmpty() )
3876 switch ( s.ShapeType() ) {
3877 case TopAbs_FACE: fmap.Add ( s ); break;
3878 case TopAbs_EDGE: emap.Add ( s ); break;
3884 if ( s.ShapeType() != TopAbs_FACE )
3885 boundarySM.push_back( sm );
3891 //================================================================================
3893 * \brief Return true if given shape is to be precomputed in order to be correctly
3894 * added to netgen mesh
3896 //================================================================================
3898 bool NETGENPlugin_Internals::isShapeToPrecompute(const TopoDS_Shape& s)
3900 int shapeID = _mesh.GetMeshDS()->ShapeToIndex( s );
3901 switch ( s.ShapeType() ) {
3902 case TopAbs_FACE : break; //return isInternalShape( shapeID ) || isBorderFace( shapeID );
3903 case TopAbs_EDGE : return isInternalEdge( shapeID );
3904 case TopAbs_VERTEX: break;
3910 //================================================================================
3912 * \brief Return SMESH
3914 //================================================================================
3916 SMESH_Mesh& NETGENPlugin_Internals::getMesh() const
3918 return const_cast<SMESH_Mesh&>( _mesh );
3921 //================================================================================
3923 * \brief Access to a counter of NETGENPlugin_NetgenLibWrapper instances
3925 //================================================================================
3927 int& NETGENPlugin_NetgenLibWrapper::instanceCounter()
3929 static int theCouner = 0;
3933 //================================================================================
3935 * \brief Initialize netgen library
3937 //================================================================================
3939 NETGENPlugin_NetgenLibWrapper::NETGENPlugin_NetgenLibWrapper()
3941 if ( instanceCounter() == 0 )
3944 ++instanceCounter();
3946 _isComputeOk = false;
3950 if ( !getenv( "KEEP_NETGEN_OUTPUT" ))
3952 // redirect all netgen output (mycout,myerr,cout) to _outputFileName
3953 _outputFileName = getOutputFileName();
3954 _ngcout = netgen::mycout;
3955 _ngcerr = netgen::myerr;
3956 netgen::mycout = new ofstream ( _outputFileName.c_str() );
3957 netgen::myerr = netgen::mycout;
3958 _coutBuffer = std::cout.rdbuf();
3960 cout << "NOTE: netgen output is redirected to file " << _outputFileName << endl;
3962 std::cout.rdbuf( netgen::mycout->rdbuf() );
3966 _ngMesh = Ng_NewMesh();
3969 //================================================================================
3971 * \brief Finish using netgen library
3973 //================================================================================
3975 NETGENPlugin_NetgenLibWrapper::~NETGENPlugin_NetgenLibWrapper()
3977 --instanceCounter();
3979 Ng_DeleteMesh( _ngMesh );
3983 std::cout.rdbuf( _coutBuffer );
3990 //================================================================================
3992 * \brief Set netgen mesh to delete at destruction
3994 //================================================================================
3996 void NETGENPlugin_NetgenLibWrapper::setMesh( Ng_Mesh* mesh )
3999 Ng_DeleteMesh( _ngMesh );
4003 //================================================================================
4005 * \brief Return a unique file name
4007 //================================================================================
4009 std::string NETGENPlugin_NetgenLibWrapper::getOutputFileName()
4011 std::string aTmpDir = SALOMEDS_Tool::GetTmpDir();
4013 TCollection_AsciiString aGenericName = (char*)aTmpDir.c_str();
4014 aGenericName += "NETGEN_";
4016 aGenericName += getpid();
4018 aGenericName += _getpid();
4020 aGenericName += "_";
4021 aGenericName += Abs((Standard_Integer)(long) aGenericName.ToCString());
4022 aGenericName += ".out";
4024 return aGenericName.ToCString();
4027 //================================================================================
4029 * \brief Remove "test.out" and "problemfaces" files in current directory
4031 //================================================================================
4033 void NETGENPlugin_NetgenLibWrapper::RemoveTmpFiles()
4035 bool rm = SMESH_File("test.out").remove() ;
4037 if ( rm && netgen::testout && instanceCounter() == 0 )
4039 delete netgen::testout;
4040 netgen::testout = 0;
4043 SMESH_File("problemfaces").remove();
4044 SMESH_File("occmesh.rep").remove();
4047 //================================================================================
4049 * \brief Remove file with netgen output
4051 //================================================================================
4053 void NETGENPlugin_NetgenLibWrapper::removeOutputFile()
4055 if ( !_outputFileName.empty() )
4059 delete netgen::mycout;
4060 netgen::mycout = _ngcout;
4061 netgen::myerr = _ngcerr;
4064 string tmpDir = SALOMEDS_Tool::GetDirFromPath ( _outputFileName );
4065 string aFileName = SALOMEDS_Tool::GetNameFromPath( _outputFileName ) + ".out";
4066 SALOMEDS::ListOfFileNames_var aFiles = new SALOMEDS::ListOfFileNames;
4068 aFiles[0] = aFileName.c_str();
4070 SALOMEDS_Tool::RemoveTemporaryFiles( tmpDir.c_str(), aFiles.in(), true );