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_LinearEdge.hxx>
36 #include <SMDS_MeshElement.hxx>
37 #include <SMDS_MeshNode.hxx>
38 #include <SMESHDS_Mesh.hxx>
39 #include <SMESH_Block.hxx>
40 #include <SMESH_Comment.hxx>
41 #include <SMESH_ComputeError.hxx>
42 #include <SMESH_File.hxx>
43 #include <SMESH_Gen_i.hxx>
44 #include <SMESH_Mesh.hxx>
45 #include <SMESH_MesherHelper.hxx>
46 #include <SMESH_subMesh.hxx>
47 #include <StdMeshers_QuadToTriaAdaptor.hxx>
48 #include <StdMeshers_ViscousLayers2D.hxx>
50 #include <SALOMEDS_Tool.hxx>
52 #include <utilities.h>
54 #include <BRepBuilderAPI_Copy.hxx>
55 #include <BRep_Tool.hxx>
56 #include <Bnd_B3d.hxx>
57 #include <NCollection_Map.hxx>
58 #include <Standard_ErrorHandler.hxx>
59 #include <Standard_ProgramError.hxx>
60 #include <TColStd_MapOfInteger.hxx>
62 #include <TopExp_Explorer.hxx>
63 #include <TopTools_DataMapIteratorOfDataMapOfShapeInteger.hxx>
64 #include <TopTools_DataMapIteratorOfDataMapOfShapeShape.hxx>
65 #include <TopTools_DataMapOfShapeInteger.hxx>
66 #include <TopTools_DataMapOfShapeShape.hxx>
67 #include <TopTools_MapOfShape.hxx>
70 // Netgen include files
74 #include <occgeom.hpp>
75 #include <meshing.hpp>
76 //#include <ngexception.hpp>
79 extern int OCCGenerateMesh (OCCGeometry&, Mesh*&, MeshingParameters&, int, int);
81 extern int OCCGenerateMesh (OCCGeometry&, Mesh*&, int, int, char*);
83 //extern void OCCSetLocalMeshSize(OCCGeometry & geom, Mesh & mesh);
84 extern MeshingParameters mparam;
85 extern volatile multithreadt multithread;
86 extern bool merge_solids;
88 // values used for occgeo.facemeshstatus
89 enum EFaceMeshStatus { FACE_NOT_TREATED = 0,
101 using namespace nglib;
105 #define nodeVec_ACCESS(index) ((SMDS_MeshNode*) nodeVec.at((index)))
107 #define nodeVec_ACCESS(index) ((SMDS_MeshNode*) nodeVec[index])
110 #define NGPOINT_COORDS(p) p(0),p(1),p(2)
113 // dump elements added to ng mesh
114 //#define DUMP_SEGMENTS
115 //#define DUMP_TRIANGLES
116 //#define DUMP_TRIANGLES_SCRIPT "/tmp/trias.py" //!< debug AddIntVerticesInSolids()
119 TopTools_IndexedMapOfShape ShapesWithLocalSize;
120 std::map<int,double> VertexId2LocalSize;
121 std::map<int,double> EdgeId2LocalSize;
122 std::map<int,double> FaceId2LocalSize;
124 //=============================================================================
128 //=============================================================================
130 NETGENPlugin_Mesher::NETGENPlugin_Mesher (SMESH_Mesh* mesh,
131 const TopoDS_Shape& aShape,
137 _fineness(NETGENPlugin_Hypothesis::GetDefaultFineness()),
138 _isViscousLayers2D(false),
147 SetDefaultParameters();
148 ShapesWithLocalSize.Clear();
149 VertexId2LocalSize.clear();
150 EdgeId2LocalSize.clear();
151 FaceId2LocalSize.clear();
154 //================================================================================
158 //================================================================================
160 NETGENPlugin_Mesher::~NETGENPlugin_Mesher()
168 //================================================================================
170 * Set pointer to NETGENPlugin_Mesher* field of the holder, that will be
171 * nullified at destruction of this
173 //================================================================================
175 void NETGENPlugin_Mesher::SetSelfPointer( NETGENPlugin_Mesher ** ptr )
186 //================================================================================
188 * \brief Initialize global NETGEN parameters with default values
190 //================================================================================
192 void NETGENPlugin_Mesher::SetDefaultParameters()
194 netgen::MeshingParameters& mparams = netgen::mparam;
195 // maximal mesh edge size
196 mparams.maxh = 0;//NETGENPlugin_Hypothesis::GetDefaultMaxSize();
198 // minimal number of segments per edge
199 mparams.segmentsperedge = NETGENPlugin_Hypothesis::GetDefaultNbSegPerEdge();
200 // rate of growth of size between elements
201 mparams.grading = NETGENPlugin_Hypothesis::GetDefaultGrowthRate();
202 // safety factor for curvatures (elements per radius)
203 mparams.curvaturesafety = NETGENPlugin_Hypothesis::GetDefaultNbSegPerRadius();
204 // create elements of second order
205 mparams.secondorder = NETGENPlugin_Hypothesis::GetDefaultSecondOrder();
206 // quad-dominated surface meshing
210 mparams.quad = NETGENPlugin_Hypothesis_2D::GetDefaultQuadAllowed();
211 _fineness = NETGENPlugin_Hypothesis::GetDefaultFineness();
212 mparams.uselocalh = NETGENPlugin_Hypothesis::GetDefaultSurfaceCurvature();
213 netgen::merge_solids = NETGENPlugin_Hypothesis::GetDefaultFuseEdges();
216 //=============================================================================
220 //=============================================================================
222 void SetLocalSize(TopoDS_Shape GeomShape, double LocalSize)
224 if ( GeomShape.IsNull() ) return;
225 TopAbs_ShapeEnum GeomType = GeomShape.ShapeType();
226 if (GeomType == TopAbs_COMPOUND) {
227 for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()) {
228 SetLocalSize(it.Value(), LocalSize);
233 if (! ShapesWithLocalSize.Contains(GeomShape))
234 key = ShapesWithLocalSize.Add(GeomShape);
236 key = ShapesWithLocalSize.FindIndex(GeomShape);
237 if (GeomType == TopAbs_VERTEX) {
238 VertexId2LocalSize[key] = LocalSize;
239 } else if (GeomType == TopAbs_EDGE) {
240 EdgeId2LocalSize[key] = LocalSize;
241 } else if (GeomType == TopAbs_FACE) {
242 FaceId2LocalSize[key] = LocalSize;
246 //=============================================================================
248 * Pass parameters to NETGEN
250 //=============================================================================
251 void NETGENPlugin_Mesher::SetParameters(const NETGENPlugin_Hypothesis* hyp)
255 netgen::MeshingParameters& mparams = netgen::mparam;
256 // Initialize global NETGEN parameters:
257 // maximal mesh segment size
258 mparams.maxh = hyp->GetMaxSize();
259 // maximal mesh element linear size
260 mparams.minh = hyp->GetMinSize();
261 // minimal number of segments per edge
262 mparams.segmentsperedge = hyp->GetNbSegPerEdge();
263 // rate of growth of size between elements
264 mparams.grading = hyp->GetGrowthRate();
265 // safety factor for curvatures (elements per radius)
266 mparams.curvaturesafety = hyp->GetNbSegPerRadius();
267 // create elements of second order
268 mparams.secondorder = hyp->GetSecondOrder() ? 1 : 0;
269 // quad-dominated surface meshing
270 mparams.quad = hyp->GetQuadAllowed() ? 1 : 0;
271 _optimize = hyp->GetOptimize();
272 _fineness = hyp->GetFineness();
273 mparams.uselocalh = hyp->GetSurfaceCurvature();
274 netgen::merge_solids = hyp->GetFuseEdges();
277 SMESH_Gen_i* smeshGen_i = SMESH_Gen_i::GetSMESHGen();
278 CORBA::Object_var anObject = smeshGen_i->GetNS()->Resolve("/myStudyManager");
279 SALOMEDS::StudyManager_var aStudyMgr = SALOMEDS::StudyManager::_narrow(anObject);
280 SALOMEDS::Study_var myStudy = aStudyMgr->GetStudyByID(hyp->GetStudyId());
282 const NETGENPlugin_Hypothesis::TLocalSize localSizes = hyp->GetLocalSizesAndEntries();
283 NETGENPlugin_Hypothesis::TLocalSize::const_iterator it = localSizes.begin();
284 for ( ; it != localSizes.end() ; it++)
286 std::string entry = (*it).first;
287 double val = (*it).second;
289 GEOM::GEOM_Object_var aGeomObj;
290 SALOMEDS::SObject_var aSObj = myStudy->FindObjectID( entry.c_str() );
291 if ( !aSObj->_is_nil() ) {
292 CORBA::Object_var obj = aSObj->GetObject();
293 aGeomObj = GEOM::GEOM_Object::_narrow(obj);
296 TopoDS_Shape S = smeshGen_i->GeomObjectToShape( aGeomObj.in() );
297 SetLocalSize(S, val);
302 //=============================================================================
304 * Pass simple parameters to NETGEN
306 //=============================================================================
308 void NETGENPlugin_Mesher::SetParameters(const NETGENPlugin_SimpleHypothesis_2D* hyp)
312 SetDefaultParameters();
315 //=============================================================================
317 * Link - a pair of integer numbers
319 //=============================================================================
323 Link(int _n1, int _n2) : n1(_n1), n2(_n2) {}
324 Link() : n1(0), n2(0) {}
325 bool Contains( int n ) const { return n == n1 || n == n2; }
326 bool IsConnected( const Link& other ) const
328 return (( Contains( other.n1 ) || Contains( other.n2 )) && ( this != &other ));
332 int HashCode(const Link& aLink, int aLimit)
334 return HashCode(aLink.n1 + aLink.n2, aLimit);
337 Standard_Boolean IsEqual(const Link& aLink1, const Link& aLink2)
339 return (( aLink1.n1 == aLink2.n1 && aLink1.n2 == aLink2.n2 ) ||
340 ( aLink1.n1 == aLink2.n2 && aLink1.n2 == aLink2.n1 ));
345 //================================================================================
347 * \brief return id of netgen point corresponding to SMDS node
349 //================================================================================
350 typedef map< const SMDS_MeshNode*, int > TNode2IdMap;
352 int ngNodeId( const SMDS_MeshNode* node,
353 netgen::Mesh& ngMesh,
354 TNode2IdMap& nodeNgIdMap)
356 int newNgId = ngMesh.GetNP() + 1;
358 TNode2IdMap::iterator node_id = nodeNgIdMap.insert( make_pair( node, newNgId )).first;
360 if ( node_id->second == newNgId)
362 #if defined(DUMP_SEGMENTS) || defined(DUMP_TRIANGLES)
363 cout << "Ng " << newNgId << " - " << node;
365 netgen::MeshPoint p( netgen::Point<3> (node->X(), node->Y(), node->Z()) );
366 ngMesh.AddPoint( p );
368 return node_id->second;
371 //================================================================================
373 * \brief Return computed EDGEs connected to the given one
375 //================================================================================
377 list< TopoDS_Edge > getConnectedEdges( const TopoDS_Edge& edge,
378 const TopoDS_Face& face,
379 const set< SMESH_subMesh* > & computedSM,
380 const SMESH_MesherHelper& helper,
381 map< SMESH_subMesh*, set< int > >& addedEdgeSM2Faces)
384 list< TopoDS_Edge > edges;
385 list< int > nbEdgesInWire;
386 /*int nbWires =*/ SMESH_Block::GetOrderedEdges( face, edges, nbEdgesInWire);
388 // find <edge> within <edges>
389 list< TopoDS_Edge >::iterator eItFwd = edges.begin();
390 for ( ; eItFwd != edges.end(); ++eItFwd )
391 if ( edge.IsSame( *eItFwd ))
393 if ( eItFwd == edges.end()) return list< TopoDS_Edge>();
395 if ( eItFwd->Orientation() >= TopAbs_INTERNAL )
397 // connected INTERNAL edges returned from GetOrderedEdges() are wrongly oriented
398 // so treat each INTERNAL edge separately
399 TopoDS_Edge e = *eItFwd;
401 edges.push_back( e );
405 // get all computed EDGEs connected to <edge>
407 list< TopoDS_Edge >::iterator eItBack = eItFwd, ePrev;
408 TopoDS_Vertex vCommon;
409 TopTools_MapOfShape eAdded; // map used not to add a seam edge twice to <edges>
412 // put edges before <edge> to <edges> back
413 while ( edges.begin() != eItFwd )
414 edges.splice( edges.end(), edges, edges.begin() );
418 while ( ++eItFwd != edges.end() )
420 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( *eItFwd );
422 bool connected = TopExp::CommonVertex( *ePrev, *eItFwd, vCommon );
423 bool computed = sm->IsMeshComputed();
424 bool added = addedEdgeSM2Faces[ sm ].count( helper.GetSubShapeID() );
425 bool doubled = !eAdded.Add( *eItFwd );
426 bool orientOK = (( ePrev ->Orientation() < TopAbs_INTERNAL ) ==
427 ( eItFwd->Orientation() < TopAbs_INTERNAL ) );
428 if ( !connected || !computed || !orientOK || added || doubled )
430 // stop advancement; move edges from tail to head
431 while ( edges.back() != *ePrev )
432 edges.splice( edges.begin(), edges, --edges.end() );
438 while ( eItBack != edges.begin() )
442 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( *eItBack );
444 bool connected = TopExp::CommonVertex( *ePrev, *eItBack, vCommon );
445 bool computed = sm->IsMeshComputed();
446 bool added = addedEdgeSM2Faces[ sm ].count( helper.GetSubShapeID() );
447 bool doubled = !eAdded.Add( *eItBack );
448 bool orientOK = (( ePrev ->Orientation() < TopAbs_INTERNAL ) ==
449 ( eItBack->Orientation() < TopAbs_INTERNAL ) );
450 if ( !connected || !computed || !orientOK || added || doubled)
453 edges.erase( edges.begin(), ePrev );
457 if ( edges.front() != edges.back() )
459 // assure that the 1st vertex is meshed
460 TopoDS_Edge eLast = edges.back();
461 while ( !SMESH_Algo::VertexNode( SMESH_MesherHelper::IthVertex( 0, edges.front()), helper.GetMeshDS())
463 edges.front() != eLast )
464 edges.splice( edges.end(), edges, edges.begin() );
469 //================================================================================
471 * \brief Make triangulation of a shape precise enough
473 //================================================================================
475 void updateTriangulation( const TopoDS_Shape& shape )
477 // static set< Poly_Triangulation* > updated;
479 // TopLoc_Location loc;
480 // TopExp_Explorer fExp( shape, TopAbs_FACE );
481 // for ( ; fExp.More(); fExp.Next() )
483 // Handle(Poly_Triangulation) triangulation =
484 // BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
485 // if ( triangulation.IsNull() ||
486 // updated.insert( triangulation.operator->() ).second )
488 // BRepTools::Clean (shape);
491 BRepMesh_IncrementalMesh e(shape, 0.01, true);
493 catch (Standard_Failure)
496 // updated.erase( triangulation.operator->() );
497 // triangulation = BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
498 // updated.insert( triangulation.operator->() );
502 //================================================================================
504 * \brief Returns a medium node either existing in SMESH of created by NETGEN
505 * \param [in] corner1 - corner node 1
506 * \param [in] corner2 - corner node 2
507 * \param [in] defaultMedium - the node created by NETGEN
508 * \param [in] helper - holder of medium nodes existing in SMESH
509 * \return const SMDS_MeshNode* - the result node
511 //================================================================================
513 const SMDS_MeshNode* mediumNode( const SMDS_MeshNode* corner1,
514 const SMDS_MeshNode* corner2,
515 const SMDS_MeshNode* defaultMedium,
516 const SMESH_MesherHelper* helper)
520 TLinkNodeMap::const_iterator l2n =
521 helper->GetTLinkNodeMap().find( SMESH_TLink( corner1, corner2 ));
522 if ( l2n != helper->GetTLinkNodeMap().end() )
523 defaultMedium = l2n->second;
525 return defaultMedium;
528 //================================================================================
530 * \brief Assure that mesh on given shapes is quadratic
532 //================================================================================
534 // void makeQuadratic( const TopTools_IndexedMapOfShape& shapes,
535 // SMESH_Mesh* mesh )
537 // for ( int i = 1; i <= shapes.Extent(); ++i )
539 // SMESHDS_SubMesh* smDS = mesh->GetMeshDS()->MeshElements( shapes(i) );
540 // if ( !smDS ) continue;
541 // SMDS_ElemIteratorPtr elemIt = smDS->GetElements();
542 // if ( !elemIt->more() ) continue;
543 // const SMDS_MeshElement* e = elemIt->next();
544 // if ( !e || e->IsQuadratic() )
547 // TIDSortedElemSet elems;
548 // elems.insert( e );
549 // while ( elemIt->more() )
550 // elems.insert( elems.end(), elemIt->next() );
552 // SMESH_MeshEditor( mesh ).ConvertToQuadratic( /*3d=*/false, elems, /*biQuad=*/false );
558 //================================================================================
560 * \brief Initialize netgen::OCCGeometry with OCCT shape
562 //================================================================================
564 void NETGENPlugin_Mesher::PrepareOCCgeometry(netgen::OCCGeometry& occgeo,
565 const TopoDS_Shape& shape,
567 list< SMESH_subMesh* > * meshedSM,
568 NETGENPlugin_Internals* intern)
570 updateTriangulation( shape );
573 BRepBndLib::Add (shape, bb);
574 double x1,y1,z1,x2,y2,z2;
575 bb.Get (x1,y1,z1,x2,y2,z2);
576 netgen::Point<3> p1 = netgen::Point<3> (x1,y1,z1);
577 netgen::Point<3> p2 = netgen::Point<3> (x2,y2,z2);
578 occgeo.boundingbox = netgen::Box<3> (p1,p2);
580 occgeo.shape = shape;
583 // fill maps of shapes of occgeo with not yet meshed subshapes
585 // get root submeshes
586 list< SMESH_subMesh* > rootSM;
587 const int shapeID = mesh.GetMeshDS()->ShapeToIndex( shape );
588 if ( shapeID > 0 ) { // SMESH_subMesh with ID 0 may exist, don't use it!
589 rootSM.push_back( mesh.GetSubMesh( shape ));
592 for ( TopoDS_Iterator it( shape ); it.More(); it.Next() )
593 rootSM.push_back( mesh.GetSubMesh( it.Value() ));
598 // add subshapes of empty submeshes
599 list< SMESH_subMesh* >::iterator rootIt = rootSM.begin(), rootEnd = rootSM.end();
600 for ( ; rootIt != rootEnd; ++rootIt ) {
601 SMESH_subMesh * root = *rootIt;
602 SMESH_subMeshIteratorPtr smIt = root->getDependsOnIterator(/*includeSelf=*/true,
603 /*complexShapeFirst=*/true);
604 // to find a right orientation of subshapes (PAL20462)
605 TopTools_IndexedMapOfShape subShapes;
606 TopExp::MapShapes(root->GetSubShape(), subShapes);
607 while ( smIt->more() )
609 SMESH_subMesh* sm = smIt->next();
610 TopoDS_Shape shape = sm->GetSubShape();
611 totNbFaces += ( shape.ShapeType() == TopAbs_FACE );
612 if ( intern && intern->isShapeToPrecompute( shape ))
614 if ( !meshedSM || sm->IsEmpty() )
616 if ( shape.ShapeType() != TopAbs_VERTEX )
617 shape = subShapes( subShapes.FindIndex( shape ));// shape -> index -> oriented shape
618 if ( shape.Orientation() >= TopAbs_INTERNAL )
619 shape.Orientation( TopAbs_FORWARD ); // isuue 0020676
620 switch ( shape.ShapeType() ) {
621 case TopAbs_FACE : occgeo.fmap.Add( shape ); break;
622 case TopAbs_EDGE : occgeo.emap.Add( shape ); break;
623 case TopAbs_VERTEX: occgeo.vmap.Add( shape ); break;
624 case TopAbs_SOLID :occgeo.somap.Add( shape ); break;
628 // collect submeshes of meshed shapes
631 const int dim = SMESH_Gen::GetShapeDim( shape );
632 meshedSM[ dim ].push_back( sm );
636 occgeo.facemeshstatus.SetSize (totNbFaces);
637 occgeo.facemeshstatus = 0;
638 occgeo.face_maxh_modified.SetSize(totNbFaces);
639 occgeo.face_maxh_modified = 0;
640 occgeo.face_maxh.SetSize(totNbFaces);
641 occgeo.face_maxh = netgen::mparam.maxh;
644 //================================================================================
646 * \brief Return a default min size value suitable for the given geometry.
648 //================================================================================
650 double NETGENPlugin_Mesher::GetDefaultMinSize(const TopoDS_Shape& geom,
651 const double maxSize)
653 updateTriangulation( geom );
657 const int* pi[4] = { &i1, &i2, &i3, &i1 };
660 TopExp_Explorer fExp( geom, TopAbs_FACE );
661 for ( ; fExp.More(); fExp.Next() )
663 Handle(Poly_Triangulation) triangulation =
664 BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
665 if ( triangulation.IsNull() ) continue;
666 const double fTol = BRep_Tool::Tolerance( TopoDS::Face( fExp.Current() ));
667 const TColgp_Array1OfPnt& points = triangulation->Nodes();
668 const Poly_Array1OfTriangle& trias = triangulation->Triangles();
669 for ( int iT = trias.Lower(); iT <= trias.Upper(); ++iT )
671 trias(iT).Get( i1, i2, i3 );
672 for ( int j = 0; j < 3; ++j )
674 double dist2 = points(*pi[j]).SquareDistance( points( *pi[j+1] ));
675 if ( dist2 < minh && fTol*fTol < dist2 )
677 bb.Add( points(*pi[j]));
681 if ( minh > 0.25 * bb.SquareExtent() ) // simple geometry, rough triangulation
683 minh = 1e-3 * sqrt( bb.SquareExtent());
684 //cout << "BND BOX minh = " <<minh << endl;
688 minh = 3 * sqrt( minh ); // triangulation for visualization is rather fine
689 //cout << "TRIANGULATION minh = " <<minh << endl;
691 if ( minh > 0.5 * maxSize )
697 //================================================================================
699 * \brief Restrict size of elements at a given point
701 //================================================================================
703 void NETGENPlugin_Mesher::RestrictLocalSize(netgen::Mesh& ngMesh,
706 const bool overrideMinH)
708 if ( size <= std::numeric_limits<double>::min() )
710 if ( netgen::mparam.minh > size )
714 ngMesh.SetMinimalH( size );
715 netgen::mparam.minh = size;
719 size = netgen::mparam.minh;
722 netgen::Point3d pi(p.X(), p.Y(), p.Z());
723 ngMesh.RestrictLocalH( pi, size );
726 //================================================================================
728 * \brief fill ngMesh with nodes and elements of computed submeshes
730 //================================================================================
732 bool NETGENPlugin_Mesher::FillNgMesh(netgen::OCCGeometry& occgeom,
733 netgen::Mesh& ngMesh,
734 vector<const SMDS_MeshNode*>& nodeVec,
735 const list< SMESH_subMesh* > & meshedSM,
736 SMESH_MesherHelper* quadHelper,
737 SMESH_ProxyMesh::Ptr proxyMesh)
739 TNode2IdMap nodeNgIdMap;
740 for ( size_t i = 1; i < nodeVec.size(); ++i )
741 nodeNgIdMap.insert( make_pair( nodeVec[i], i ));
743 TopTools_MapOfShape visitedShapes;
744 map< SMESH_subMesh*, set< int > > visitedEdgeSM2Faces;
745 set< SMESH_subMesh* > computedSM( meshedSM.begin(), meshedSM.end() );
747 SMESH_MesherHelper helper (*_mesh);
749 int faceNgID = ngMesh.GetNFD();
751 list< SMESH_subMesh* >::const_iterator smIt, smEnd = meshedSM.end();
752 for ( smIt = meshedSM.begin(); smIt != smEnd; ++smIt )
754 SMESH_subMesh* sm = *smIt;
755 if ( !visitedShapes.Add( sm->GetSubShape() ))
758 const SMESHDS_SubMesh * smDS = sm->GetSubMeshDS();
759 if ( !smDS ) continue;
761 switch ( sm->GetSubShape().ShapeType() )
763 case TopAbs_EDGE: { // EDGE
764 // ----------------------
765 TopoDS_Edge geomEdge = TopoDS::Edge( sm->GetSubShape() );
766 if ( geomEdge.Orientation() >= TopAbs_INTERNAL )
767 geomEdge.Orientation( TopAbs_FORWARD ); // issue 0020676
769 // Add ng segments for each not meshed FACE the EDGE bounds
770 PShapeIteratorPtr fIt = helper.GetAncestors( geomEdge, *sm->GetFather(), TopAbs_FACE );
771 while ( const TopoDS_Shape * anc = fIt->next() )
773 faceNgID = occgeom.fmap.FindIndex( *anc );
775 continue; // meshed face
777 int faceSMDSId = helper.GetMeshDS()->ShapeToIndex( *anc );
778 if ( visitedEdgeSM2Faces[ sm ].count( faceSMDSId ))
779 continue; // already treated EDGE
781 TopoDS_Face face = TopoDS::Face( occgeom.fmap( faceNgID ));
782 if ( face.Orientation() >= TopAbs_INTERNAL )
783 face.Orientation( TopAbs_FORWARD ); // issue 0020676
785 // get all meshed EDGEs of the FACE connected to geomEdge (issue 0021140)
786 helper.SetSubShape( face );
787 list< TopoDS_Edge > edges = getConnectedEdges( geomEdge, face, computedSM, helper,
788 visitedEdgeSM2Faces );
790 continue; // wrong ancestor?
792 // find out orientation of <edges> within <face>
793 TopoDS_Edge eNotSeam = edges.front();
794 if ( helper.HasSeam() )
796 list< TopoDS_Edge >::iterator eIt = edges.begin();
797 while ( helper.IsRealSeam( *eIt )) ++eIt;
798 if ( eIt != edges.end() )
801 TopAbs_Orientation fOri = helper.GetSubShapeOri( face, eNotSeam );
802 bool isForwad = ( fOri == eNotSeam.Orientation() || fOri >= TopAbs_INTERNAL );
804 // get all nodes from connected <edges>
805 const bool isQuad = smDS->IsQuadratic();
806 StdMeshers_FaceSide fSide( face, edges, _mesh, isForwad, isQuad );
807 const vector<UVPtStruct>& points = fSide.GetUVPtStruct();
808 if ( points.empty() )
809 return false; // invalid node params?
810 int i, nbSeg = fSide.NbSegments();
812 // remember EDGEs of fSide to treat only once
813 for ( int iE = 0; iE < fSide.NbEdges(); ++iE )
814 visitedEdgeSM2Faces[ helper.GetMesh()->GetSubMesh( fSide.Edge(iE )) ].insert(faceSMDSId);
816 double otherSeamParam = 0;
821 int prevNgId = ngNodeId( points[0].node, ngMesh, nodeNgIdMap );
823 for ( i = 0; i < nbSeg; ++i )
825 const UVPtStruct& p1 = points[ i ];
826 const UVPtStruct& p2 = points[ i+1 ];
828 if ( p1.node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX ) //an EDGE begins
831 if ( helper.IsRealSeam( p1.node->getshapeId() ))
833 TopoDS_Edge e = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
834 isSeam = helper.IsRealSeam( e );
837 otherSeamParam = helper.GetOtherParam( helper.GetPeriodicIndex() & 1 ? p2.u : p2.v );
844 seg[1] = prevNgId = ngNodeId( p2.node, ngMesh, nodeNgIdMap );
845 // node param on curve
846 seg.epgeominfo[ 0 ].dist = p1.param;
847 seg.epgeominfo[ 1 ].dist = p2.param;
849 seg.epgeominfo[ 0 ].u = p1.u;
850 seg.epgeominfo[ 0 ].v = p1.v;
851 seg.epgeominfo[ 1 ].u = p2.u;
852 seg.epgeominfo[ 1 ].v = p2.v;
854 //geomEdge = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
855 //seg.epgeominfo[ 0 ].edgenr = seg.epgeominfo[ 1 ].edgenr = occgeom.emap.FindIndex( geomEdge );
857 //seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
858 seg.si = faceNgID; // = geom.fmap.FindIndex (face);
859 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
860 ngMesh.AddSegment (seg);
862 SMESH_TNodeXYZ np1( p1.node ), np2( p2.node );
863 RestrictLocalSize( ngMesh, 0.5*(np1+np2), (np1-np2).Modulus() );
866 cout << "Segment: " << seg.edgenr << " on SMESH face " << helper.GetMeshDS()->ShapeToIndex( face ) << endl
867 << "\tface index: " << seg.si << endl
868 << "\tp1: " << seg[0] << endl
869 << "\tp2: " << seg[1] << endl
870 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
871 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
872 //<< "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
873 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
874 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl;
875 //<< "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
879 if ( helper.GetPeriodicIndex() && 1 ) {
880 seg.epgeominfo[ 0 ].u = otherSeamParam;
881 seg.epgeominfo[ 1 ].u = otherSeamParam;
882 swap (seg.epgeominfo[0].v, seg.epgeominfo[1].v);
884 seg.epgeominfo[ 0 ].v = otherSeamParam;
885 seg.epgeominfo[ 1 ].v = otherSeamParam;
886 swap (seg.epgeominfo[0].u, seg.epgeominfo[1].u);
888 swap( seg[0], seg[1] );
889 swap( seg.epgeominfo[0].dist, seg.epgeominfo[1].dist );
890 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
891 ngMesh.AddSegment( seg );
893 cout << "Segment: " << seg.edgenr << endl
894 << "\t is SEAM (reverse) of the previous. "
895 << " Other " << (helper.GetPeriodicIndex() && 1 ? "U" : "V")
896 << " = " << otherSeamParam << endl;
899 else if ( fOri == TopAbs_INTERNAL )
901 swap( seg[0], seg[1] );
902 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
903 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
904 ngMesh.AddSegment( seg );
906 cout << "Segment: " << seg.edgenr << endl << "\t is REVERSE of the previous" << endl;
910 } // loop on geomEdge ancestors
912 if ( quadHelper ) // remember medium nodes of sub-meshes
914 SMDS_ElemIteratorPtr edges = smDS->GetElements();
915 while ( edges->more() )
917 const SMDS_MeshElement* e = edges->next();
918 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshEdge*>( e )))
924 } // case TopAbs_EDGE
926 case TopAbs_FACE: { // FACE
927 // ----------------------
928 const TopoDS_Face& geomFace = TopoDS::Face( sm->GetSubShape() );
929 helper.SetSubShape( geomFace );
930 bool isInternalFace = ( geomFace.Orientation() == TopAbs_INTERNAL );
932 // Find solids the geomFace bounds
933 int solidID1 = 0, solidID2 = 0;
934 StdMeshers_QuadToTriaAdaptor* quadAdaptor =
935 dynamic_cast<StdMeshers_QuadToTriaAdaptor*>( proxyMesh.get() );
938 solidID1 = occgeom.somap.FindIndex( quadAdaptor->GetShape() );
942 PShapeIteratorPtr solidIt = helper.GetAncestors( geomFace, *sm->GetFather(), TopAbs_SOLID);
943 while ( const TopoDS_Shape * solid = solidIt->next() )
945 int id = occgeom.somap.FindIndex ( *solid );
946 if ( solidID1 && id != solidID1 ) solidID2 = id;
950 // Add ng face descriptors of meshed faces
952 ngMesh.AddFaceDescriptor( netgen::FaceDescriptor( faceNgID, solidID1, solidID2, 0 ));
954 // if second oreder is required, even already meshed faces must be passed to NETGEN
955 int fID = occgeom.fmap.Add( geomFace );
956 if ( occgeom.facemeshstatus.Size() < fID ) occgeom.facemeshstatus.SetSize( fID );
957 occgeom.facemeshstatus[ fID-1 ] = netgen::FACE_MESHED_OK;
958 while ( fID < faceNgID ) // geomFace is already in occgeom.fmap, add a copy
960 fID = occgeom.fmap.Add( BRepBuilderAPI_Copy( geomFace, /*copyGeom=*/false ));
961 if ( occgeom.facemeshstatus.Size() < fID ) occgeom.facemeshstatus.SetSize( fID );
962 occgeom.facemeshstatus[ fID-1 ] = netgen::FACE_MESHED_OK;
964 // Problem with the second order in a quadrangular mesh remains.
965 // 1) All quadrangles generated by NETGEN are moved to an inexistent face
966 // by FillSMesh() (find "AddFaceDescriptor")
967 // 2) Temporary triangles generated by StdMeshers_QuadToTriaAdaptor
968 // are on faces where quadrangles were.
969 // Due to these 2 points, wrong geom faces are used while conversion to quadratic
970 // of the mentioned above quadrangles and triangles
972 // Orient the face correctly in solidID1 (issue 0020206)
973 bool reverse = false;
975 TopoDS_Shape solid = occgeom.somap( solidID1 );
976 TopAbs_Orientation faceOriInSolid = helper.GetSubShapeOri( solid, geomFace );
977 if ( faceOriInSolid >= 0 )
979 helper.IsReversedSubMesh( TopoDS::Face( geomFace.Oriented( faceOriInSolid )));
982 // Add surface elements
984 netgen::Element2d tri(3);
985 tri.SetIndex( faceNgID );
986 SMESH_TNodeXYZ xyz[3];
988 #ifdef DUMP_TRIANGLES
989 cout << "SMESH face " << helper.GetMeshDS()->ShapeToIndex( geomFace )
990 << " internal="<<isInternalFace << endl;
993 smDS = proxyMesh->GetSubMesh( geomFace );
995 SMDS_ElemIteratorPtr faces = smDS->GetElements();
996 while ( faces->more() )
998 const SMDS_MeshElement* f = faces->next();
999 if ( f->NbNodes() % 3 != 0 ) // not triangle
1001 PShapeIteratorPtr solidIt=helper.GetAncestors(geomFace,*sm->GetFather(),TopAbs_SOLID);
1002 if ( const TopoDS_Shape * solid = solidIt->next() )
1003 sm = _mesh->GetSubMesh( *solid );
1004 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1005 smError.reset( new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH,"Not triangle sub-mesh"));
1006 smError->myBadElements.push_back( f );
1010 for ( int i = 0; i < 3; ++i )
1012 const SMDS_MeshNode* node = f->GetNode( i ), * inFaceNode=0;
1015 // get node UV on face
1016 int shapeID = node->getshapeId();
1017 if ( helper.IsSeamShape( shapeID ))
1019 if ( helper.IsSeamShape( f->GetNodeWrap( i+1 )->getshapeId() ))
1020 inFaceNode = f->GetNodeWrap( i-1 );
1022 inFaceNode = f->GetNodeWrap( i+1 );
1024 gp_XY uv = helper.GetNodeUV( geomFace, node, inFaceNode );
1026 int ind = reverse ? 3-i : i+1;
1027 tri.GeomInfoPi(ind).u = uv.X();
1028 tri.GeomInfoPi(ind).v = uv.Y();
1029 tri.PNum (ind) = ngNodeId( node, ngMesh, nodeNgIdMap );
1032 // pass a triangle size to NG size-map
1033 double size = ( ( xyz[0] - xyz[1] ).Modulus() +
1034 ( xyz[1] - xyz[2] ).Modulus() +
1035 ( xyz[2] - xyz[0] ).Modulus() ) / 3;
1036 gp_XYZ gc = ( xyz[0] + xyz[1] + xyz[2] ) / 3;
1037 RestrictLocalSize( ngMesh, gc, size, /*overrideMinH=*/false );
1039 ngMesh.AddSurfaceElement (tri);
1040 #ifdef DUMP_TRIANGLES
1041 cout << tri << endl;
1044 if ( isInternalFace )
1046 swap( tri[1], tri[2] );
1047 ngMesh.AddSurfaceElement (tri);
1048 #ifdef DUMP_TRIANGLES
1049 cout << tri << endl;
1054 if ( quadHelper ) // remember medium nodes of sub-meshes
1056 SMDS_ElemIteratorPtr faces = smDS->GetElements();
1057 while ( faces->more() )
1059 const SMDS_MeshElement* f = faces->next();
1060 if ( !quadHelper->AddTLinks( static_cast< const SMDS_MeshFace*>( f )))
1066 } // case TopAbs_FACE
1068 case TopAbs_VERTEX: { // VERTEX
1069 // --------------------------
1070 // issue 0021405. Add node only if a VERTEX is shared by a not meshed EDGE,
1071 // else netgen removes a free node and nodeVector becomes invalid
1072 PShapeIteratorPtr ansIt = helper.GetAncestors( sm->GetSubShape(),
1076 while ( const TopoDS_Shape* e = ansIt->next() )
1078 SMESH_subMesh* eSub = helper.GetMesh()->GetSubMesh( *e );
1079 if (( toAdd = ( eSub->IsEmpty() && !SMESH_Algo::isDegenerated( TopoDS::Edge( *e )))))
1084 SMDS_NodeIteratorPtr nodeIt = smDS->GetNodes();
1085 if ( nodeIt->more() )
1086 ngNodeId( nodeIt->next(), ngMesh, nodeNgIdMap );
1092 } // loop on submeshes
1095 nodeVec.resize( ngMesh.GetNP() + 1 );
1096 TNode2IdMap::iterator node_NgId, nodeNgIdEnd = nodeNgIdMap.end();
1097 for ( node_NgId = nodeNgIdMap.begin(); node_NgId != nodeNgIdEnd; ++node_NgId)
1098 nodeVec[ node_NgId->second ] = node_NgId->first;
1103 //================================================================================
1105 * \brief Duplicate mesh faces on internal geom faces
1107 //================================================================================
1109 void NETGENPlugin_Mesher::FixIntFaces(const netgen::OCCGeometry& occgeom,
1110 netgen::Mesh& ngMesh,
1111 NETGENPlugin_Internals& internalShapes)
1113 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1115 // find ng indices of internal faces
1117 for ( int ngFaceID = 1; ngFaceID <= occgeom.fmap.Extent(); ++ngFaceID )
1119 int smeshID = meshDS->ShapeToIndex( occgeom.fmap( ngFaceID ));
1120 if ( internalShapes.isInternalShape( smeshID ))
1121 ngFaceIds.insert( ngFaceID );
1123 if ( !ngFaceIds.empty() )
1126 int i, nbFaces = ngMesh.GetNSE();
1127 for ( i = 1; i <= nbFaces; ++i)
1129 netgen::Element2d elem = ngMesh.SurfaceElement(i);
1130 if ( ngFaceIds.count( elem.GetIndex() ))
1132 swap( elem[1], elem[2] );
1133 ngMesh.AddSurfaceElement (elem);
1139 //================================================================================
1141 * \brief Tries to heal the mesh on a FACE. The FACE is supposed to be partially
1142 * meshed due to NETGEN failure
1143 * \param [in] occgeom - geometry
1144 * \param [in,out] ngMesh - the mesh to fix
1145 * \param [inout] faceID - ID of the FACE to fix the mesh on
1146 * \return bool - is mesh is or becomes OK
1148 //================================================================================
1150 bool NETGENPlugin_Mesher::FixFaceMesh(const netgen::OCCGeometry& occgeom,
1151 netgen::Mesh& ngMesh,
1154 // we address a case where the FACE is almost fully meshed except small holes
1155 // of usually triangular shape at FACE boundary (IPAL52861)
1157 // The case appeared to be not simple: holes only look triangular but
1158 // indeed are a self intersecting polygon. A reason of the bug was in coincident
1159 // NG points on a seam edge. But the code below is very nice, leave it for
1164 if ( occgeom.fmap.Extent() < faceID )
1166 //const TopoDS_Face& face = TopoDS::Face( occgeom.fmap( faceID ));
1168 // find free links on the FACE
1169 NCollection_Map<Link> linkMap;
1170 for ( int iF = 1; iF <= ngMesh.GetNSE(); ++iF )
1172 const netgen::Element2d& elem = ngMesh.SurfaceElement(iF);
1173 if ( faceID != elem.GetIndex() )
1175 int n0 = elem[ elem.GetNP() - 1 ];
1176 for ( int i = 0; i < elem.GetNP(); ++i )
1179 Link link( n0, n1 );
1180 if ( !linkMap.Add( link ))
1181 linkMap.Remove( link );
1185 // add/remove boundary links
1186 for ( int iSeg = 1; iSeg <= ngMesh.GetNSeg(); ++iSeg )
1188 const netgen::Segment& seg = ngMesh.LineSegment( iSeg );
1189 if ( seg.si != faceID ) // !edgeIDs.Contains( seg.edgenr ))
1191 Link link( seg[1], seg[0] ); // reverse!!!
1192 if ( !linkMap.Add( link ))
1193 linkMap.Remove( link );
1195 if ( linkMap.IsEmpty() )
1197 if ( linkMap.Extent() < 3 )
1200 // make triangles of the links
1202 netgen::Element2d tri(3);
1203 tri.SetIndex ( faceID );
1205 NCollection_Map<Link>::Iterator linkIt( linkMap );
1206 Link link1 = linkIt.Value();
1207 // look for a link connected to link1
1208 NCollection_Map<Link>::Iterator linkIt2 = linkIt;
1209 for ( linkIt2.Next(); linkIt2.More(); linkIt2.Next() )
1211 const Link& link2 = linkIt2.Value();
1212 if ( link2.IsConnected( link1 ))
1214 // look for a link connected to both link1 and link2
1215 NCollection_Map<Link>::Iterator linkIt3 = linkIt2;
1216 for ( linkIt3.Next(); linkIt3.More(); linkIt3.Next() )
1218 const Link& link3 = linkIt3.Value();
1219 if ( link3.IsConnected( link1 ) &&
1220 link3.IsConnected( link2 ) )
1225 tri[2] = ( link2.Contains( link1.n1 ) ? link2.n1 : link3.n1 );
1226 if ( tri[0] == tri[2] || tri[1] == tri[2] )
1228 ngMesh.AddSurfaceElement( tri );
1230 // prepare for the next tria search
1231 if ( linkMap.Extent() == 3 )
1233 linkMap.Remove( link3 );
1234 linkMap.Remove( link2 );
1236 linkMap.Remove( link1 );
1237 link1 = linkIt.Value();
1250 //================================================================================
1251 // define gp_XY_Subtracted pointer to function calling gp_XY::Subtracted(gp_XY)
1252 gp_XY_FunPtr(Subtracted);
1253 //gp_XY_FunPtr(Added);
1255 //================================================================================
1257 * \brief Evaluate distance between two 2d points along the surface
1259 //================================================================================
1261 double evalDist( const gp_XY& uv1,
1263 const Handle(Geom_Surface)& surf,
1264 const int stopHandler=-1)
1266 if ( stopHandler > 0 ) // continue recursion
1268 gp_XY mid = SMESH_MesherHelper::GetMiddleUV( surf, uv1, uv2 );
1269 return evalDist( uv1,mid, surf, stopHandler-1 ) + evalDist( mid,uv2, surf, stopHandler-1 );
1271 double dist3D = surf->Value( uv1.X(), uv1.Y() ).Distance( surf->Value( uv2.X(), uv2.Y() ));
1272 if ( stopHandler == 0 ) // stop recursion
1275 // start recursion if necessary
1276 double dist2D = SMESH_MesherHelper::ApplyIn2D(surf, uv1, uv2, gp_XY_Subtracted, 0).Modulus();
1277 if ( fabs( dist3D - dist2D ) < dist2D * 1e-10 )
1278 return dist3D; // equal parametrization of a planar surface
1280 return evalDist( uv1, uv2, surf, 3 ); // start recursion
1283 //================================================================================
1285 * \brief Data of vertex internal in geom face
1287 //================================================================================
1291 gp_XY uv; //!< UV in face parametric space
1292 int ngId; //!< ng id of corrsponding node
1293 gp_XY uvClose; //!< UV of closest boundary node
1294 int ngIdClose; //!< ng id of closest boundary node
1297 //================================================================================
1299 * \brief Data of vertex internal in solid
1301 //================================================================================
1305 int ngId; //!< ng id of corresponding node
1306 int ngIdClose; //!< ng id of closest 2d mesh element
1307 int ngIdCloseN; //!< ng id of closest node of the closest 2d mesh element
1310 inline double dist2(const netgen::MeshPoint& p1, const netgen::MeshPoint& p2)
1312 return gp_Pnt( NGPOINT_COORDS(p1)).SquareDistance( gp_Pnt( NGPOINT_COORDS(p2)));
1316 //================================================================================
1318 * \brief Make netgen take internal vertices in faces into account by adding
1319 * segments including internal vertices
1321 * This function works in supposition that 1D mesh is already computed in ngMesh
1323 //================================================================================
1325 void NETGENPlugin_Mesher::AddIntVerticesInFaces(const netgen::OCCGeometry& occgeom,
1326 netgen::Mesh& ngMesh,
1327 vector<const SMDS_MeshNode*>& nodeVec,
1328 NETGENPlugin_Internals& internalShapes)
1330 if ((int) nodeVec.size() < ngMesh.GetNP() )
1331 nodeVec.resize( ngMesh.GetNP(), 0 );
1333 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1334 SMESH_MesherHelper helper( internalShapes.getMesh() );
1336 const map<int,list<int> >& face2Vert = internalShapes.getFacesWithVertices();
1337 map<int,list<int> >::const_iterator f2v = face2Vert.begin();
1338 for ( ; f2v != face2Vert.end(); ++f2v )
1340 const TopoDS_Face& face = TopoDS::Face( meshDS->IndexToShape( f2v->first ));
1341 if ( face.IsNull() ) continue;
1342 int faceNgID = occgeom.fmap.FindIndex (face);
1343 if ( faceNgID < 0 ) continue;
1345 TopLoc_Location loc;
1346 Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc);
1348 helper.SetSubShape( face );
1349 helper.SetElementsOnShape( true );
1351 // Get data of internal vertices and add them to ngMesh
1353 multimap< double, TIntVData > dist2VData; // sort vertices by distance from boundary nodes
1355 int i, nbSegInit = ngMesh.GetNSeg();
1357 // boundary characteristics
1358 double totSegLen2D = 0;
1361 const list<int>& iVertices = f2v->second;
1362 list<int>::const_iterator iv = iVertices.begin();
1363 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1366 // get node on vertex
1367 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1368 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1371 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1372 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1373 nV = SMESH_Algo::VertexNode( V, meshDS );
1374 if ( !nV ) continue;
1377 netgen::MeshPoint mp( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1378 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1379 vData.ngId = ngMesh.GetNP();
1380 nodeVec.push_back( nV );
1384 vData.uv = helper.GetNodeUV( face, nV, 0, &uvOK );
1385 if ( !uvOK ) helper.CheckNodeUV( face, nV, vData.uv, BRep_Tool::Tolerance(V),/*force=*/1);
1387 // loop on all segments of the face to find the node closest to vertex and to count
1388 // average segment 2d length
1389 double closeDist2 = numeric_limits<double>::max(), dist2;
1391 for (i = 1; i <= ngMesh.GetNSeg(); ++i)
1393 netgen::Segment & seg = ngMesh.LineSegment(i);
1394 if ( seg.si != faceNgID ) continue;
1396 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1398 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1399 if ( ngIdLast == seg[ iEnd ] ) continue;
1400 dist2 = helper.ApplyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1401 if ( dist2 < closeDist2 )
1402 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1403 ngIdLast = seg[ iEnd ];
1407 totSegLen2D += helper.ApplyIn2D(surf, uv[0], uv[1], gp_XY_Subtracted, false).Modulus();
1411 dist2VData.insert( make_pair( closeDist2, vData ));
1414 if ( totNbSeg == 0 ) break;
1415 double avgSegLen2d = totSegLen2D / totNbSeg;
1417 // Loop on vertices to add segments
1419 multimap< double, TIntVData >::iterator dist_vData = dist2VData.begin();
1420 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1422 double closeDist2 = dist_vData->first, dist2;
1423 TIntVData & vData = dist_vData->second;
1425 // try to find more close node among segments added for internal vertices
1426 for (i = nbSegInit+1; i <= ngMesh.GetNSeg(); ++i)
1428 netgen::Segment & seg = ngMesh.LineSegment(i);
1429 if ( seg.si != faceNgID ) continue;
1431 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1433 uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1434 dist2 = helper.ApplyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1435 if ( dist2 < closeDist2 )
1436 vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1439 // decide whether to use the closest node as the second end of segment or to
1440 // create a new point
1441 int segEnd1 = vData.ngId;
1442 int segEnd2 = vData.ngIdClose; // to use closest node
1443 gp_XY uvV = vData.uv, uvP = vData.uvClose;
1444 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1445 double nodeDist2D = sqrt( closeDist2 );
1446 double nodeDist3D = evalDist( vData.uv, vData.uvClose, surf );
1447 bool avgLenOK = ( avgSegLen2d < 0.75 * nodeDist2D );
1448 bool hintLenOK = ( segLenHint < 0.75 * nodeDist3D );
1449 //cout << "uvV " << uvV.X() <<","<<uvV.Y() << " ";
1450 if ( hintLenOK || avgLenOK )
1452 // create a point between the closest node and V
1455 double r = min( 0.5, ( hintLenOK ? segLenHint/nodeDist3D : avgSegLen2d/nodeDist2D ));
1456 // direction from V to closet node in 2D
1457 gp_Dir2d v2n( helper.ApplyIn2D(surf, uvP, uvV, gp_XY_Subtracted, false ));
1459 uvP = vData.uv + r * nodeDist2D * v2n.XY();
1460 gp_Pnt P = surf->Value( uvP.X(), uvP.Y() ).Transformed( loc );
1462 netgen::MeshPoint mp( netgen::Point<3> (P.X(), P.Y(), P.Z()));
1463 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1464 segEnd2 = ngMesh.GetNP();
1465 //cout << "Middle " << r << " uv " << uvP.X() << "," << uvP.Y() << "( " << ngMesh.Point(segEnd2).X()<<","<<ngMesh.Point(segEnd2).Y()<<","<<ngMesh.Point(segEnd2).Z()<<" )"<< endl;
1466 SMDS_MeshNode * nP = helper.AddNode(P.X(), P.Y(), P.Z());
1467 nodeVec.push_back( nP );
1469 //else cout << "at Node " << " uv " << uvP.X() << "," << uvP.Y() << endl;
1472 netgen::Segment seg;
1474 if ( segEnd1 > segEnd2 ) swap( segEnd1, segEnd2 ), swap( uvV, uvP );
1475 seg[0] = segEnd1; // ng node id
1476 seg[1] = segEnd2; // ng node id
1477 seg.edgenr = ngMesh.GetNSeg() + 1;// segment id
1480 seg.epgeominfo[ 0 ].dist = 0; // param on curve
1481 seg.epgeominfo[ 0 ].u = uvV.X();
1482 seg.epgeominfo[ 0 ].v = uvV.Y();
1483 seg.epgeominfo[ 1 ].dist = 1; // param on curve
1484 seg.epgeominfo[ 1 ].u = uvP.X();
1485 seg.epgeominfo[ 1 ].v = uvP.Y();
1487 // seg.epgeominfo[ 0 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1488 // seg.epgeominfo[ 1 ].edgenr = 10; // = geom.emap.FindIndex(edge);
1490 ngMesh.AddSegment (seg);
1492 // add reverse segment
1493 swap( seg[0], seg[1] );
1494 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1495 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1496 ngMesh.AddSegment (seg);
1502 //================================================================================
1504 * \brief Make netgen take internal vertices in solids into account by adding
1505 * faces including internal vertices
1507 * This function works in supposition that 2D mesh is already computed in ngMesh
1509 //================================================================================
1511 void NETGENPlugin_Mesher::AddIntVerticesInSolids(const netgen::OCCGeometry& occgeom,
1512 netgen::Mesh& ngMesh,
1513 vector<const SMDS_MeshNode*>& nodeVec,
1514 NETGENPlugin_Internals& internalShapes)
1516 #ifdef DUMP_TRIANGLES_SCRIPT
1517 // create a python script making a mesh containing triangles added for internal vertices
1518 ofstream py(DUMP_TRIANGLES_SCRIPT);
1519 py << "import SMESH"<< endl
1520 << "from salome.smesh import smeshBuilder"<<endl
1521 << "smesh = smeshBuilder.New(salome.myStudy)"<<endl
1522 << "m = smesh.Mesh(name='triangles')" << endl;
1524 if ((int) nodeVec.size() < ngMesh.GetNP() )
1525 nodeVec.resize( ngMesh.GetNP(), 0 );
1527 SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1528 SMESH_MesherHelper helper( internalShapes.getMesh() );
1530 const map<int,list<int> >& so2Vert = internalShapes.getSolidsWithVertices();
1531 map<int,list<int> >::const_iterator s2v = so2Vert.begin();
1532 for ( ; s2v != so2Vert.end(); ++s2v )
1534 const TopoDS_Shape& solid = meshDS->IndexToShape( s2v->first );
1535 if ( solid.IsNull() ) continue;
1536 int solidNgID = occgeom.somap.FindIndex (solid);
1537 if ( solidNgID < 0 && !occgeom.somap.IsEmpty() ) continue;
1539 helper.SetSubShape( solid );
1540 helper.SetElementsOnShape( true );
1542 // find ng indices of faces within the solid
1544 for (TopExp_Explorer fExp(solid, TopAbs_FACE); fExp.More(); fExp.Next() )
1545 ngFaceIds.insert( occgeom.fmap.FindIndex( fExp.Current() ));
1546 if ( ngFaceIds.size() == 1 && *ngFaceIds.begin() == 0 )
1547 ngFaceIds.insert( 1 );
1549 // Get data of internal vertices and add them to ngMesh
1551 multimap< double, TIntVSoData > dist2VData; // sort vertices by distance from ng faces
1553 int i, nbFaceInit = ngMesh.GetNSE();
1555 // boundary characteristics
1556 double totSegLen = 0;
1559 const list<int>& iVertices = s2v->second;
1560 list<int>::const_iterator iv = iVertices.begin();
1561 for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1564 const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1566 // get node on vertex
1567 const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1570 SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1571 sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1572 nV = SMESH_Algo::VertexNode( V, meshDS );
1573 if ( !nV ) continue;
1576 netgen::MeshPoint mpV( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1577 ngMesh.AddPoint ( mpV, 1, netgen::FIXEDPOINT );
1578 vData.ngId = ngMesh.GetNP();
1579 nodeVec.push_back( nV );
1581 // loop on all 2d elements to find the one closest to vertex and to count
1582 // average segment length
1583 double closeDist2 = numeric_limits<double>::max(), avgDist2;
1584 for (i = 1; i <= ngMesh.GetNSE(); ++i)
1586 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1587 if ( !ngFaceIds.count( elem.GetIndex() )) continue;
1589 multimap< double, int> dist2nID; // sort nodes of element by distance from V
1590 for ( int j = 0; j < elem.GetNP(); ++j)
1592 netgen::MeshPoint mp = ngMesh.Point( elem[j] );
1593 double d2 = dist2( mpV, mp );
1594 dist2nID.insert( make_pair( d2, elem[j] ));
1595 avgDist2 += d2 / elem.GetNP();
1597 totNbSeg++, totSegLen+= sqrt( dist2( mp, ngMesh.Point( elem[(j+1)%elem.GetNP()])));
1599 double dist = dist2nID.begin()->first; //avgDist2;
1600 if ( dist < closeDist2 )
1601 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= dist;
1603 dist2VData.insert( make_pair( closeDist2, vData ));
1606 if ( totNbSeg == 0 ) break;
1607 double avgSegLen = totSegLen / totNbSeg;
1609 // Loop on vertices to add triangles
1611 multimap< double, TIntVSoData >::iterator dist_vData = dist2VData.begin();
1612 for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1614 double closeDist2 = dist_vData->first;
1615 TIntVSoData & vData = dist_vData->second;
1617 const netgen::MeshPoint& mpV = ngMesh.Point( vData.ngId );
1619 // try to find more close face among ones added for internal vertices
1620 for (i = nbFaceInit+1; i <= ngMesh.GetNSE(); ++i)
1622 double avgDist2 = 0;
1623 multimap< double, int> dist2nID;
1624 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1625 for ( int j = 0; j < elem.GetNP(); ++j)
1627 double d = dist2( mpV, ngMesh.Point( elem[j] ));
1628 dist2nID.insert( make_pair( d, elem[j] ));
1629 avgDist2 += d / elem.GetNP();
1630 if ( avgDist2 < closeDist2 )
1631 vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= avgDist2;
1634 // sort nodes of the closest face by angle with vector from V to the closest node
1635 const double tol = numeric_limits<double>::min();
1636 map< double, int > angle2ID;
1637 const netgen::Element2d& closeFace = ngMesh.SurfaceElement( vData.ngIdClose );
1638 netgen::MeshPoint mp[2];
1639 mp[0] = ngMesh.Point( vData.ngIdCloseN );
1640 gp_XYZ p1( NGPOINT_COORDS( mp[0] ));
1641 gp_XYZ pV( NGPOINT_COORDS( mpV ));
1642 gp_Vec v2p1( pV, p1 );
1643 double distN1 = v2p1.Magnitude();
1644 if ( distN1 <= tol ) continue;
1646 for ( int j = 0; j < closeFace.GetNP(); ++j)
1648 mp[1] = ngMesh.Point( closeFace[j] );
1649 gp_Vec v2p( pV, gp_Pnt( NGPOINT_COORDS( mp[1] )) );
1650 angle2ID.insert( make_pair( v2p1.Angle( v2p ), closeFace[j]));
1652 // get node with angle of 60 degrees or greater
1653 map< double, int >::iterator angle_id = angle2ID.lower_bound( 60. * M_PI / 180. );
1654 if ( angle_id == angle2ID.end() ) angle_id = --angle2ID.end();
1655 const double minAngle = 30. * M_PI / 180.;
1656 const double angle = angle_id->first;
1657 bool angleOK = ( angle > minAngle );
1659 // find points to create a triangle
1660 netgen::Element2d tri(3);
1662 tri[0] = vData.ngId;
1663 tri[1] = vData.ngIdCloseN; // to use the closest nodes
1664 tri[2] = angle_id->second; // to use the node with best angle
1666 // decide whether to use the closest node and the node with best angle or to create new ones
1667 for ( int isBestAngleN = 0; isBestAngleN < 2; ++isBestAngleN )
1669 bool createNew = !angleOK; //, distOK = true;
1671 int triInd = isBestAngleN ? 2 : 1;
1672 mp[isBestAngleN] = ngMesh.Point( tri[triInd] );
1677 double distN2 = sqrt( dist2( mpV, mp[isBestAngleN]));
1678 createNew = ( fabs( distN2 - distN1 ) > 0.25 * distN1 );
1680 else if ( angle < tol )
1682 v2p1.SetX( v2p1.X() + 1e-3 );
1688 double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1689 bool avgLenOK = ( avgSegLen < 0.75 * distN1 );
1690 bool hintLenOK = ( segLenHint < 0.75 * distN1 );
1691 createNew = (createNew || avgLenOK || hintLenOK );
1692 // we create a new node not closer than 0.5 to the closest face
1693 // in order not to clash with other close face
1694 double r = min( 0.5, ( hintLenOK ? segLenHint : avgSegLen ) / distN1 );
1695 distFromV = r * distN1;
1699 // create a new point, between the node and the vertex if angleOK
1700 gp_XYZ p( NGPOINT_COORDS( mp[isBestAngleN] ));
1701 gp_Vec v2p( pV, p ); v2p.Normalize();
1702 if ( isBestAngleN && !angleOK )
1703 p = p1 + gp_Dir( v2p.XYZ() - v2p1.XYZ()).XYZ() * distN1 * 0.95;
1705 p = pV + v2p.XYZ() * distFromV;
1707 if ( !isBestAngleN ) p1 = p, distN1 = distFromV;
1709 mp[isBestAngleN].SetPoint( netgen::Point<3> (p.X(), p.Y(), p.Z()));
1710 ngMesh.AddPoint ( mp[isBestAngleN], 1, netgen::SURFACEPOINT );
1711 tri[triInd] = ngMesh.GetNP();
1712 nodeVec.push_back( helper.AddNode( p.X(), p.Y(), p.Z()) );
1715 ngMesh.AddSurfaceElement (tri);
1716 swap( tri[1], tri[2] );
1717 ngMesh.AddSurfaceElement (tri);
1719 #ifdef DUMP_TRIANGLES_SCRIPT
1720 py << "n1 = m.AddNode( "<< mpV(0)<<", "<< mpV(1)<<", "<< mpV(2)<<") "<< endl
1721 << "n2 = m.AddNode( "<< mp[0](0)<<", "<< mp[0](1)<<", "<< mp[0](2)<<") "<< endl
1722 << "n3 = m.AddNode( "<< mp[1](0)<<", "<< mp[1](1)<<", "<< mp[1](2)<<" )" << endl
1723 << "m.AddFace([n1,n2,n3])" << endl;
1725 } // loop on internal vertices of a solid
1727 } // loop on solids with internal vertices
1730 //================================================================================
1732 * \brief Fill netgen mesh with segments of a FACE
1733 * \param ngMesh - netgen mesh
1734 * \param geom - container of OCCT geometry to mesh
1735 * \param wires - data of nodes on FACE boundary
1736 * \param helper - mesher helper holding the FACE
1737 * \param nodeVec - vector of nodes in which node index == netgen ID
1738 * \retval SMESH_ComputeErrorPtr - error description
1740 //================================================================================
1742 SMESH_ComputeErrorPtr
1743 NETGENPlugin_Mesher::AddSegmentsToMesh(netgen::Mesh& ngMesh,
1744 netgen::OCCGeometry& geom,
1745 const TSideVector& wires,
1746 SMESH_MesherHelper& helper,
1747 vector< const SMDS_MeshNode* > & nodeVec,
1748 const bool overrideMinH)
1750 // ----------------------------
1751 // Check wires and count nodes
1752 // ----------------------------
1754 for ( size_t iW = 0; iW < wires.size(); ++iW )
1756 StdMeshers_FaceSidePtr wire = wires[ iW ];
1757 if ( wire->MissVertexNode() )
1759 // Commented for issue 0020960. It worked for the case, let's wait for case where it doesn't.
1760 // It seems that there is no reason for this limitation
1762 // (new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH, "Missing nodes on vertices"));
1764 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1765 if ((int) uvPtVec.size() != wire->NbPoints() )
1766 return SMESH_ComputeError::New(COMPERR_BAD_INPUT_MESH,
1767 SMESH_Comment("Unexpected nb of points on wire ") << iW
1768 << ": " << uvPtVec.size()<<" != "<<wire->NbPoints());
1769 nbNodes += wire->NbPoints();
1771 nodeVec.reserve( nodeVec.size() + nbNodes + 1 );
1772 if ( nodeVec.empty() )
1773 nodeVec.push_back( 0 );
1775 // -----------------
1777 // -----------------
1779 const bool wasNgMeshEmpty = ( ngMesh.GetNP() < 1 ); /* true => this method is called by
1780 NETGENPlugin_NETGEN_2D_ONLY */
1782 // map for nodes on vertices since they can be shared between wires
1783 // ( issue 0020676, face_int_box.brep) and nodes built by NETGEN
1784 map<const SMDS_MeshNode*, int > node2ngID;
1785 if ( !wasNgMeshEmpty ) // fill node2ngID with nodes built by NETGEN
1787 set< int > subIDs; // ids of sub-shapes of the FACE
1788 for ( size_t iW = 0; iW < wires.size(); ++iW )
1790 StdMeshers_FaceSidePtr wire = wires[ iW ];
1791 for ( int iE = 0, nbE = wire->NbEdges(); iE < nbE; ++iE )
1793 subIDs.insert( wire->EdgeID( iE ));
1794 subIDs.insert( helper.GetMeshDS()->ShapeToIndex( wire->FirstVertex( iE )));
1797 for ( size_t ngID = 1; ngID < nodeVec.size(); ++ngID )
1798 if ( subIDs.count( nodeVec[ngID]->getshapeId() ))
1799 node2ngID.insert( make_pair( nodeVec[ngID], ngID ));
1802 const int solidID = 0, faceID = geom.fmap.FindIndex( helper.GetSubShape() );
1803 if ( ngMesh.GetNFD() < 1 )
1804 ngMesh.AddFaceDescriptor( netgen::FaceDescriptor( faceID, solidID, solidID, 0 ));
1806 for ( size_t iW = 0; iW < wires.size(); ++iW )
1808 StdMeshers_FaceSidePtr wire = wires[ iW ];
1809 const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1810 const int nbSegments = wire->NbPoints() - 1;
1812 // assure the 1st node to be in node2ngID, which is needed to correctly
1813 // "close chain of segments" (see below) in case if the 1st node is not
1814 // onVertex because it is on a Viscous layer
1815 node2ngID.insert( make_pair( uvPtVec[ 0 ].node, ngMesh.GetNP() + 1 ));
1817 // compute length of every segment
1818 vector<double> segLen( nbSegments );
1819 for ( int i = 0; i < nbSegments; ++i )
1820 segLen[i] = SMESH_TNodeXYZ( uvPtVec[ i ].node ).Distance( uvPtVec[ i+1 ].node );
1822 int edgeID = 1, posID = -2;
1823 bool isInternalWire = false;
1824 double vertexNormPar = 0;
1825 const int prevNbNGSeg = ngMesh.GetNSeg();
1826 for ( int i = 0; i < nbSegments; ++i ) // loop on segments
1828 // Add the first point of a segment
1830 const SMDS_MeshNode * n = uvPtVec[ i ].node;
1831 const int posShapeID = n->getshapeId();
1832 bool onVertex = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX );
1833 bool onEdge = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_EDGE );
1835 // skip nodes on degenerated edges
1836 if ( helper.IsDegenShape( posShapeID ) &&
1837 helper.IsDegenShape( uvPtVec[ i+1 ].node->getshapeId() ))
1840 int ngID1 = ngMesh.GetNP() + 1, ngID2 = ngID1+1;
1841 if ( onVertex || ( !wasNgMeshEmpty && onEdge ) || helper.IsRealSeam( posShapeID ))
1842 ngID1 = node2ngID.insert( make_pair( n, ngID1 )).first->second;
1843 if ( ngID1 > ngMesh.GetNP() )
1845 netgen::MeshPoint mp( netgen::Point<3> (n->X(), n->Y(), n->Z()) );
1846 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1847 nodeVec.push_back( n );
1849 else // n is in ngMesh already, and ngID2 in prev segment is wrong
1851 ngID2 = ngMesh.GetNP() + 1;
1852 if ( i > 0 ) // prev segment belongs to same wire
1854 netgen::Segment& prevSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
1861 netgen::Segment seg;
1863 seg[0] = ngID1; // ng node id
1864 seg[1] = ngID2; // ng node id
1865 seg.edgenr = ngMesh.GetNSeg() + 1; // ng segment id
1866 seg.si = faceID; // = geom.fmap.FindIndex (face);
1868 for ( int iEnd = 0; iEnd < 2; ++iEnd)
1870 const UVPtStruct& pnt = uvPtVec[ i + iEnd ];
1872 seg.epgeominfo[ iEnd ].dist = pnt.param; // param on curve
1873 seg.epgeominfo[ iEnd ].u = pnt.u;
1874 seg.epgeominfo[ iEnd ].v = pnt.v;
1876 // find out edge id and node parameter on edge
1877 onVertex = ( pnt.normParam + 1e-10 > vertexNormPar );
1878 if ( onVertex || posShapeID != posID )
1881 double normParam = pnt.normParam;
1883 normParam = 0.5 * ( uvPtVec[ i ].normParam + uvPtVec[ i+1 ].normParam );
1884 int edgeIndexInWire = wire->EdgeIndex( normParam );
1885 vertexNormPar = wire->LastParameter( edgeIndexInWire );
1886 const TopoDS_Edge& edge = wire->Edge( edgeIndexInWire );
1887 edgeID = geom.emap.FindIndex( edge );
1889 isInternalWire = ( edge.Orientation() == TopAbs_INTERNAL );
1890 // if ( onVertex ) // param on curve is different on each of two edges
1891 // seg.epgeominfo[ iEnd ].dist = helper.GetNodeU( edge, pnt.node );
1893 seg.epgeominfo[ iEnd ].edgenr = edgeID; // = geom.emap.FindIndex(edge);
1896 ngMesh.AddSegment (seg);
1898 // restrict size of elements near the segment
1899 SMESH_TNodeXYZ np1( n ), np2( uvPtVec[ i+1 ].node );
1900 // get an average size of adjacent segments to avoid sharp change of
1901 // element size (regression on issue 0020452, note 0010898)
1902 int iPrev = SMESH_MesherHelper::WrapIndex( i-1, nbSegments );
1903 int iNext = SMESH_MesherHelper::WrapIndex( i+1, nbSegments );
1904 double sumH = segLen[ iPrev ] + segLen[ i ] + segLen[ iNext ];
1905 int nbSeg = ( int( segLen[ iPrev ] > sumH / 100.) +
1906 int( segLen[ i ] > sumH / 100.) +
1907 int( segLen[ iNext ] > sumH / 100.));
1909 RestrictLocalSize( ngMesh, 0.5*(np1+np2), sumH / nbSeg, overrideMinH );
1911 if ( isInternalWire )
1913 swap (seg[0], seg[1]);
1914 swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1915 seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1916 ngMesh.AddSegment (seg);
1918 } // loop on segments on a wire
1920 // close chain of segments
1921 if ( nbSegments > 0 )
1923 netgen::Segment& lastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() - int( isInternalWire ));
1924 const SMDS_MeshNode * lastNode = uvPtVec.back().node;
1925 lastSeg[1] = node2ngID.insert( make_pair( lastNode, lastSeg[1] )).first->second;
1926 if ( lastSeg[1] > ngMesh.GetNP() )
1928 netgen::MeshPoint mp( netgen::Point<3> (lastNode->X(), lastNode->Y(), lastNode->Z()) );
1929 ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1930 nodeVec.push_back( lastNode );
1932 if ( isInternalWire )
1934 netgen::Segment& realLastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
1935 realLastSeg[0] = lastSeg[1];
1939 #ifdef DUMP_SEGMENTS
1940 cout << "BEGIN WIRE " << iW << endl;
1941 for ( int i = prevNbNGSeg+1; i <= ngMesh.GetNSeg(); ++i )
1943 netgen::Segment& seg = ngMesh.LineSegment( i );
1945 netgen::Segment& prevSeg = ngMesh.LineSegment( i-1 );
1946 if ( seg[0] == prevSeg[1] && seg[1] == prevSeg[0] )
1948 cout << "Segment: " << seg.edgenr << endl << "\tis REVRESE of the previous one" << endl;
1952 cout << "Segment: " << seg.edgenr << endl
1953 << "\tp1: " << seg[0] << " n" << nodeVec[ seg[0]]->GetID() << endl
1954 << "\tp2: " << seg[1] << " n" << nodeVec[ seg[1]]->GetID() << endl
1955 << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
1956 << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
1957 << "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
1958 << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
1959 << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl
1960 << "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
1962 cout << "--END WIRE " << iW << endl;
1964 SMESH_Comment __not_unused_variable( prevNbNGSeg );
1967 } // loop on WIREs of a FACE
1969 // add a segment instead of an internal vertex
1970 if ( wasNgMeshEmpty )
1972 NETGENPlugin_Internals intShapes( *helper.GetMesh(), helper.GetSubShape(), /*is3D=*/false );
1973 AddIntVerticesInFaces( geom, ngMesh, nodeVec, intShapes );
1975 ngMesh.CalcSurfacesOfNode();
1980 //================================================================================
1982 * \brief Fill SMESH mesh according to contents of netgen mesh
1983 * \param occgeo - container of OCCT geometry to mesh
1984 * \param ngMesh - netgen mesh
1985 * \param initState - bn of entities in netgen mesh before computing
1986 * \param sMesh - SMESH mesh to fill in
1987 * \param nodeVec - vector of nodes in which node index == netgen ID
1988 * \param comment - returns problem description
1989 * \param quadHelper - holder of medium nodes of sub-meshes
1990 * \retval int - error
1992 //================================================================================
1994 int NETGENPlugin_Mesher::FillSMesh(const netgen::OCCGeometry& occgeo,
1995 netgen::Mesh& ngMesh,
1996 const NETGENPlugin_ngMeshInfo& initState,
1998 std::vector<const SMDS_MeshNode*>& nodeVec,
1999 SMESH_Comment& comment,
2000 SMESH_MesherHelper* quadHelper)
2002 int nbNod = ngMesh.GetNP();
2003 int nbSeg = ngMesh.GetNSeg();
2004 int nbFac = ngMesh.GetNSE();
2005 int nbVol = ngMesh.GetNE();
2007 SMESHDS_Mesh* meshDS = sMesh.GetMeshDS();
2009 // quadHelper is used for either
2010 // 1) making quadratic elements when a lower dimention mesh is loaded
2011 // to SMESH before convertion to quadratic by NETGEN
2012 // 2) sewing of quadratic elements with quadratic elements of sub-meshes
2013 if ( quadHelper && !quadHelper->GetIsQuadratic() && quadHelper->GetTLinkNodeMap().empty() )
2016 // -------------------------------------
2017 // Create and insert nodes into nodeVec
2018 // -------------------------------------
2020 nodeVec.resize( nbNod + 1 );
2021 int i, nbInitNod = initState._nbNodes;
2022 for (i = nbInitNod+1; i <= nbNod; ++i )
2024 const netgen::MeshPoint& ngPoint = ngMesh.Point(i);
2025 SMDS_MeshNode* node = NULL;
2026 TopoDS_Vertex aVert;
2027 // First, netgen creates nodes on vertices in occgeo.vmap,
2028 // so node index corresponds to vertex index
2029 // but (issue 0020776) netgen does not create nodes with equal coordinates
2030 if ( i-nbInitNod <= occgeo.vmap.Extent() )
2032 gp_Pnt p ( NGPOINT_COORDS(ngPoint) );
2033 for (int iV = i-nbInitNod; aVert.IsNull() && iV <= occgeo.vmap.Extent(); ++iV)
2035 aVert = TopoDS::Vertex( occgeo.vmap( iV ));
2036 gp_Pnt pV = BRep_Tool::Pnt( aVert );
2037 if ( p.SquareDistance( pV ) > 1e-20 )
2040 node = const_cast<SMDS_MeshNode*>( SMESH_Algo::VertexNode( aVert, meshDS ));
2043 if (!node) // node not found on vertex
2045 node = meshDS->AddNode( NGPOINT_COORDS( ngPoint ));
2046 if (!aVert.IsNull())
2047 meshDS->SetNodeOnVertex(node, aVert);
2052 // -------------------------------------------
2053 // Create mesh segments along geometric edges
2054 // -------------------------------------------
2056 int nbInitSeg = initState._nbSegments;
2057 for (i = nbInitSeg+1; i <= nbSeg; ++i )
2059 const netgen::Segment& seg = ngMesh.LineSegment(i);
2061 int pinds[3] = { seg.pnums[0], seg.pnums[1], seg.pnums[2] };
2064 for (int j=0; j < 3; ++j)
2066 int pind = pinds[j];
2067 if (pind <= 0 || !nodeVec_ACCESS(pind))
2075 int aGeomEdgeInd = seg.epgeominfo[j].edgenr;
2076 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
2077 aEdge = TopoDS::Edge(occgeo.emap(aGeomEdgeInd));
2079 param = seg.epgeominfo[j].dist;
2082 else // middle point
2084 param = param2 * 0.5;
2086 if (!aEdge.IsNull() && nodeVec_ACCESS(pind)->getshapeId() < 1)
2088 meshDS->SetNodeOnEdge(nodeVec_ACCESS(pind), aEdge, param);
2093 SMDS_MeshEdge* edge = 0;
2094 if (nbp == 2) // second order ?
2096 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1])))
2098 if ( quadHelper ) // final mesh must be quadratic
2099 edge = quadHelper->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
2101 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
2105 if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
2106 nodeVec_ACCESS(pinds[2])))
2108 edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
2109 nodeVec_ACCESS(pinds[2]));
2113 if ( comment.empty() ) comment << "Cannot create a mesh edge";
2114 MESSAGE("Cannot create a mesh edge");
2115 nbSeg = nbFac = nbVol = 0;
2118 if ( !aEdge.IsNull() && edge->getshapeId() < 1 )
2119 meshDS->SetMeshElementOnShape(edge, aEdge);
2121 else if ( comment.empty() )
2123 comment << "Invalid netgen segment #" << i;
2127 // ----------------------------------------
2128 // Create mesh faces along geometric faces
2129 // ----------------------------------------
2131 int nbInitFac = initState._nbFaces;
2132 int quadFaceID = ngMesh.GetNFD() + 1;
2133 if ( nbInitFac < nbFac )
2134 // add a faces descriptor to exclude qudrangle elements generated by NETGEN
2135 // from computation of 3D mesh
2136 ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(quadFaceID, /*solid1=*/0, /*solid2=*/0, 0));
2138 vector<const SMDS_MeshNode*> nodes;
2139 for (i = nbInitFac+1; i <= nbFac; ++i )
2141 const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
2142 const int aGeomFaceInd = elem.GetIndex();
2144 if (aGeomFaceInd > 0 && aGeomFaceInd <= occgeo.fmap.Extent())
2145 aFace = TopoDS::Face(occgeo.fmap(aGeomFaceInd));
2147 for ( int j = 1; j <= elem.GetNP(); ++j )
2149 int pind = elem.PNum(j);
2150 if ( pind < 1 || pind >= (int) nodeVec.size() )
2152 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind))
2154 nodes.push_back( node );
2155 if (!aFace.IsNull() && node->getshapeId() < 1)
2157 const netgen::PointGeomInfo& pgi = elem.GeomInfoPi(j);
2158 meshDS->SetNodeOnFace(node, aFace, pgi.u, pgi.v);
2162 if ((int) nodes.size() != elem.GetNP() )
2164 if ( comment.empty() )
2165 comment << "Invalid netgen 2d element #" << i;
2166 continue; // bad node ids
2168 SMDS_MeshFace* face = NULL;
2169 switch (elem.GetType())
2172 if ( quadHelper ) // final mesh must be quadratic
2173 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2]);
2175 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2]);
2178 if ( quadHelper ) // final mesh must be quadratic
2179 face = quadHelper->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2181 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
2182 // exclude qudrangle elements from computation of 3D mesh
2183 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2186 nodes[5] = mediumNode( nodes[0],nodes[1],nodes[5], quadHelper );
2187 nodes[3] = mediumNode( nodes[1],nodes[2],nodes[3], quadHelper );
2188 nodes[4] = mediumNode( nodes[2],nodes[0],nodes[4], quadHelper );
2189 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[5],nodes[3],nodes[4]);
2192 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2193 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2194 nodes[5] = mediumNode( nodes[2],nodes[3],nodes[5], quadHelper );
2195 nodes[6] = mediumNode( nodes[3],nodes[0],nodes[6], quadHelper );
2196 face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3],
2197 nodes[4],nodes[7],nodes[5],nodes[6]);
2198 // exclude qudrangle elements from computation of 3D mesh
2199 const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
2202 MESSAGE("NETGEN created a face of unexpected type, ignoring");
2207 if ( comment.empty() ) comment << "Cannot create a mesh face";
2208 MESSAGE("Cannot create a mesh face");
2209 nbSeg = nbFac = nbVol = 0;
2212 if ( !aFace.IsNull() )
2213 meshDS->SetMeshElementOnShape( face, aFace );
2216 // ------------------
2217 // Create tetrahedra
2218 // ------------------
2220 for ( i = 1; i <= nbVol; ++i )
2222 const netgen::Element& elem = ngMesh.VolumeElement(i);
2223 int aSolidInd = elem.GetIndex();
2224 TopoDS_Solid aSolid;
2225 if ( aSolidInd > 0 && aSolidInd <= occgeo.somap.Extent() )
2226 aSolid = TopoDS::Solid(occgeo.somap(aSolidInd));
2228 for ( int j = 1; j <= elem.GetNP(); ++j )
2230 int pind = elem.PNum(j);
2231 if ( pind < 1 || pind >= (int)nodeVec.size() )
2233 if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind) )
2235 nodes.push_back(node);
2236 if ( !aSolid.IsNull() && node->getshapeId() < 1 )
2237 meshDS->SetNodeInVolume(node, aSolid);
2240 if ((int) nodes.size() != elem.GetNP() )
2242 if ( comment.empty() )
2243 comment << "Invalid netgen 3d element #" << i;
2246 SMDS_MeshVolume* vol = NULL;
2247 switch ( elem.GetType() )
2250 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3]);
2253 nodes[4] = mediumNode( nodes[0],nodes[1],nodes[4], quadHelper );
2254 nodes[7] = mediumNode( nodes[1],nodes[2],nodes[7], quadHelper );
2255 nodes[5] = mediumNode( nodes[2],nodes[0],nodes[5], quadHelper );
2256 nodes[6] = mediumNode( nodes[0],nodes[3],nodes[6], quadHelper );
2257 nodes[8] = mediumNode( nodes[1],nodes[3],nodes[8], quadHelper );
2258 nodes[9] = mediumNode( nodes[2],nodes[3],nodes[9], quadHelper );
2259 vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3],
2260 nodes[4],nodes[7],nodes[5],nodes[6],nodes[8],nodes[9]);
2263 MESSAGE("NETGEN created a volume of unexpected type, ignoring");
2268 if ( comment.empty() ) comment << "Cannot create a mesh volume";
2269 MESSAGE("Cannot create a mesh volume");
2270 nbSeg = nbFac = nbVol = 0;
2273 if (!aSolid.IsNull())
2274 meshDS->SetMeshElementOnShape(vol, aSolid);
2276 return comment.empty() ? 0 : 1;
2281 //================================================================================
2283 * \brief Restrict size of elements on the given edge
2285 //================================================================================
2287 void setLocalSize(const TopoDS_Edge& edge,
2291 if ( size <= std::numeric_limits<double>::min() )
2293 Standard_Real u1, u2;
2294 Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, u1, u2);
2295 if ( curve.IsNull() )
2297 TopoDS_Iterator vIt( edge );
2298 if ( !vIt.More() ) return;
2299 gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex( vIt.Value() ));
2300 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
2304 const int nb = (int)( 1.5 * SMESH_Algo::EdgeLength( edge ) / size );
2305 Standard_Real delta = (u2-u1)/nb;
2306 for(int i=0; i<nb; i++)
2308 Standard_Real u = u1 + delta*i;
2309 gp_Pnt p = curve->Value(u);
2310 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
2311 netgen::Point3d pi(p.X(), p.Y(), p.Z());
2312 double resultSize = mesh.GetH(pi);
2313 if ( resultSize - size > 0.1*size )
2314 // netgen does restriction iff oldH/newH > 1.2 (localh.cpp:136)
2315 NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), resultSize/1.201 );
2320 //================================================================================
2322 * \brief Convert error into text
2324 //================================================================================
2326 std::string text(int err)
2331 SMESH_Comment("Error in netgen::OCCGenerateMesh() at ") << netgen::multithread.task;
2334 //================================================================================
2336 * \brief Convert exception into text
2338 //================================================================================
2340 std::string text(Standard_Failure& ex)
2342 SMESH_Comment str("Exception in netgen::OCCGenerateMesh()");
2343 str << " at " << netgen::multithread.task
2344 << ": " << ex.DynamicType()->Name();
2345 if ( ex.GetMessageString() && strlen( ex.GetMessageString() ))
2346 str << ": " << ex.GetMessageString();
2349 //================================================================================
2351 * \brief Convert exception into text
2353 //================================================================================
2355 std::string text(netgen::NgException& ex)
2357 SMESH_Comment str("NgException");
2358 if ( strlen( netgen::multithread.task ) > 0 )
2359 str << " at " << netgen::multithread.task;
2360 str << ": " << ex.What();
2364 //================================================================================
2366 * \brief Looks for triangles lying on a SOLID
2368 //================================================================================
2370 bool hasBadElemOnSolid( const list<const SMDS_MeshElement*>& elems,
2371 SMESH_subMesh* solidSM )
2373 TopTools_IndexedMapOfShape solidSubs;
2374 TopExp::MapShapes( solidSM->GetSubShape(), solidSubs );
2375 SMESHDS_Mesh* mesh = solidSM->GetFather()->GetMeshDS();
2377 list<const SMDS_MeshElement*>::const_iterator e = elems.begin();
2378 for ( ; e != elems.end(); ++e )
2380 const SMDS_MeshElement* elem = *e;
2381 // if ( elem->GetType() != SMDSAbs_Face ) -- 23047
2383 int nbNodesOnSolid = 0, nbNodes = elem->NbNodes();
2384 SMDS_NodeIteratorPtr nIt = elem->nodeIterator();
2385 while ( nIt->more() )
2387 const SMDS_MeshNode* n = nIt->next();
2388 const TopoDS_Shape& s = mesh->IndexToShape( n->getshapeId() );
2389 nbNodesOnSolid += ( !s.IsNull() && solidSubs.Contains( s ));
2390 if ( nbNodesOnSolid > 2 ||
2391 nbNodesOnSolid == nbNodes)
2398 const double edgeMeshingTime = 0.001;
2399 const double faceMeshingTime = 0.019;
2400 const double edgeFaceMeshingTime = edgeMeshingTime + faceMeshingTime;
2401 const double faceOptimizTime = 0.06;
2402 const double voluMeshingTime = 0.15;
2403 const double volOptimizeTime = 0.77;
2406 //=============================================================================
2408 * Here we are going to use the NETGEN mesher
2410 //=============================================================================
2412 bool NETGENPlugin_Mesher::Compute()
2414 NETGENPlugin_NetgenLibWrapper ngLib;
2416 netgen::MeshingParameters& mparams = netgen::mparam;
2417 MESSAGE("Compute with:\n"
2418 " max size = " << mparams.maxh << "\n"
2419 " segments per edge = " << mparams.segmentsperedge);
2421 " growth rate = " << mparams.grading << "\n"
2422 " elements per radius = " << mparams.curvaturesafety << "\n"
2423 " second order = " << mparams.secondorder << "\n"
2424 " quad allowed = " << mparams.quad << "\n"
2425 " surface curvature = " << mparams.uselocalh << "\n"
2426 " fuse edges = " << netgen::merge_solids);
2428 SMESH_ComputeErrorPtr error = SMESH_ComputeError::New();
2429 SMESH_MesherHelper quadHelper( *_mesh );
2430 quadHelper.SetIsQuadratic( mparams.secondorder );
2432 static string debugFile = "/tmp/ngMesh.py"; /* to call toPython( _ngMesh, debugFile )
2433 while debugging netgen */
2434 // -------------------------
2435 // Prepare OCC geometry
2436 // -------------------------
2438 netgen::OCCGeometry occgeo;
2439 list< SMESH_subMesh* > meshedSM[3]; // for 0-2 dimensions
2440 NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
2441 PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
2444 _totalTime = edgeFaceMeshingTime;
2446 _totalTime += faceOptimizTime;
2448 _totalTime += voluMeshingTime + ( _optimize ? volOptimizeTime : 0 );
2449 double doneTime = 0;
2452 _curShapeIndex = -1;
2454 // -------------------------
2455 // Generate the mesh
2456 // -------------------------
2459 NETGENPlugin_ngMeshInfo initState; // it remembers size of ng mesh equal to size of Smesh
2461 SMESH_Comment comment;
2464 // vector of nodes in which node index == netgen ID
2465 vector< const SMDS_MeshNode* > nodeVec;
2473 // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
2474 mparams.uselocalh = false;
2475 mparams.grading = 0.8; // not limitited size growth
2477 if ( _simpleHyp->GetNumberOfSegments() )
2479 mparams.maxh = occgeo.boundingbox.Diam();
2482 mparams.maxh = _simpleHyp->GetLocalLength();
2485 if ( mparams.maxh == 0.0 )
2486 mparams.maxh = occgeo.boundingbox.Diam();
2487 if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
2488 mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
2490 // Local size on faces
2491 occgeo.face_maxh = mparams.maxh;
2493 // Let netgen create _ngMesh and calculate element size on not meshed shapes
2497 int startWith = netgen::MESHCONST_ANALYSE;
2498 int endWith = netgen::MESHCONST_ANALYSE;
2503 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2505 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2507 if(netgen::multithread.terminate)
2510 comment << text(err);
2512 catch (Standard_Failure& ex)
2514 comment << text(ex);
2516 err = 0; //- MESHCONST_ANALYSE isn't so important step
2519 ngLib.setMesh(( Ng_Mesh*) _ngMesh );
2521 _ngMesh->ClearFaceDescriptors(); // we make descriptors our-self
2525 // Pass 1D simple parameters to NETGEN
2526 // --------------------------------
2527 int nbSeg = _simpleHyp->GetNumberOfSegments();
2528 double segSize = _simpleHyp->GetLocalLength();
2529 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
2531 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
2533 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
2534 setLocalSize( e, segSize, *_ngMesh );
2537 else // if ( ! _simpleHyp )
2539 // Local size on vertices and edges
2540 // --------------------------------
2541 for(std::map<int,double>::const_iterator it=EdgeId2LocalSize.begin(); it!=EdgeId2LocalSize.end(); it++)
2543 int key = (*it).first;
2544 double hi = (*it).second;
2545 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2546 const TopoDS_Edge& e = TopoDS::Edge(shape);
2547 setLocalSize( e, hi, *_ngMesh );
2549 for(std::map<int,double>::const_iterator it=VertexId2LocalSize.begin(); it!=VertexId2LocalSize.end(); it++)
2551 int key = (*it).first;
2552 double hi = (*it).second;
2553 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2554 const TopoDS_Vertex& v = TopoDS::Vertex(shape);
2555 gp_Pnt p = BRep_Tool::Pnt(v);
2556 NETGENPlugin_Mesher::RestrictLocalSize( *_ngMesh, p.XYZ(), hi );
2558 for(map<int,double>::const_iterator it=FaceId2LocalSize.begin();
2559 it!=FaceId2LocalSize.end(); it++)
2561 int key = (*it).first;
2562 double val = (*it).second;
2563 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2564 int faceNgID = occgeo.fmap.FindIndex(shape);
2565 occgeo.SetFaceMaxH(faceNgID, val);
2566 for ( TopExp_Explorer edgeExp( shape, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
2567 setLocalSize( TopoDS::Edge( edgeExp.Current() ), val, *_ngMesh );
2571 // Precompute internal edges (issue 0020676) in order to
2572 // add mesh on them correctly (twice) to netgen mesh
2573 if ( !err && internals.hasInternalEdges() )
2575 // load internal shapes into OCCGeometry
2576 netgen::OCCGeometry intOccgeo;
2577 internals.getInternalEdges( intOccgeo.fmap, intOccgeo.emap, intOccgeo.vmap, meshedSM );
2578 intOccgeo.boundingbox = occgeo.boundingbox;
2579 intOccgeo.shape = occgeo.shape;
2580 intOccgeo.face_maxh.SetSize(intOccgeo.fmap.Extent());
2581 intOccgeo.face_maxh = netgen::mparam.maxh;
2582 netgen::Mesh *tmpNgMesh = NULL;
2586 // compute local H on internal shapes in the main mesh
2587 //OCCSetLocalMeshSize(intOccgeo, *_ngMesh); it deletes _ngMesh->localH
2589 // let netgen create a temporary mesh
2591 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2593 netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2595 if(netgen::multithread.terminate)
2598 // copy LocalH from the main to temporary mesh
2599 initState.transferLocalH( _ngMesh, tmpNgMesh );
2601 // compute mesh on internal edges
2602 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2604 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2606 err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2608 comment << text(err);
2610 catch (Standard_Failure& ex)
2612 comment << text(ex);
2615 initState.restoreLocalH( tmpNgMesh );
2617 // fill SMESH by netgen mesh
2618 vector< const SMDS_MeshNode* > tmpNodeVec;
2619 FillSMesh( intOccgeo, *tmpNgMesh, initState, *_mesh, tmpNodeVec, comment );
2620 err = ( err || !comment.empty() );
2622 nglib::Ng_DeleteMesh((nglib::Ng_Mesh*)tmpNgMesh);
2625 // Fill _ngMesh with nodes and segments of computed submeshes
2628 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_0D ]) &&
2629 FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_1D ], &quadHelper));
2631 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2636 startWith = endWith = netgen::MESHCONST_MESHEDGES;
2641 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2643 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2645 if(netgen::multithread.terminate)
2648 comment << text(err);
2650 catch (Standard_Failure& ex)
2652 comment << text(ex);
2657 _ticTime = ( doneTime += edgeMeshingTime ) / _totalTime / _progressTic;
2659 mparams.uselocalh = true; // restore as it is used at surface optimization
2661 // ---------------------
2662 // compute surface mesh
2663 // ---------------------
2666 // Pass 2D simple parameters to NETGEN
2668 if ( double area = _simpleHyp->GetMaxElementArea() ) {
2670 mparams.maxh = sqrt(2. * area/sqrt(3.0));
2671 mparams.grading = 0.4; // moderate size growth
2674 // length from edges
2675 if ( _ngMesh->GetNSeg() ) {
2676 double edgeLength = 0;
2677 TopTools_MapOfShape visitedEdges;
2678 for ( TopExp_Explorer exp( _shape, TopAbs_EDGE ); exp.More(); exp.Next() )
2679 if( visitedEdges.Add(exp.Current()) )
2680 edgeLength += SMESH_Algo::EdgeLength( TopoDS::Edge( exp.Current() ));
2681 // we have to multiply length by 2 since for each TopoDS_Edge there
2682 // are double set of NETGEN edges, in other words, we have to
2683 // divide _ngMesh->GetNSeg() by 2.
2684 mparams.maxh = 2*edgeLength / _ngMesh->GetNSeg();
2687 mparams.maxh = 1000;
2689 mparams.grading = 0.2; // slow size growth
2691 mparams.quad = _simpleHyp->GetAllowQuadrangles();
2692 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2693 _ngMesh->SetGlobalH (mparams.maxh);
2694 netgen::Box<3> bb = occgeo.GetBoundingBox();
2695 bb.Increase (bb.Diam()/20);
2696 _ngMesh->SetLocalH (bb.PMin(), bb.PMax(), mparams.grading);
2699 // Care of vertices internal in faces (issue 0020676)
2700 if ( internals.hasInternalVertexInFace() )
2702 // store computed segments in SMESH in order not to create SMESH
2703 // edges for ng segments added by AddIntVerticesInFaces()
2704 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2705 // add segments to faces with internal vertices
2706 AddIntVerticesInFaces( occgeo, *_ngMesh, nodeVec, internals );
2707 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2710 // Build viscous layers
2711 if ( _isViscousLayers2D ||
2712 StdMeshers_ViscousLayers2D::HasProxyMesh( TopoDS::Face( occgeo.fmap(1) ), *_mesh ))
2714 if ( !internals.hasInternalVertexInFace() ) {
2715 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2716 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2718 SMESH_ProxyMesh::Ptr viscousMesh;
2719 SMESH_MesherHelper helper( *_mesh );
2720 for ( int faceID = 1; faceID <= occgeo.fmap.Extent(); ++faceID )
2722 const TopoDS_Face& F = TopoDS::Face( occgeo.fmap( faceID ));
2723 viscousMesh = StdMeshers_ViscousLayers2D::Compute( *_mesh, F );
2726 if ( viscousMesh->NbProxySubMeshes() == 0 )
2728 // exclude from computation ng segments built on EDGEs of F
2729 for (int i = 1; i <= _ngMesh->GetNSeg(); i++)
2731 netgen::Segment & seg = _ngMesh->LineSegment(i);
2732 if (seg.si == faceID)
2735 // add new segments to _ngMesh instead of excluded ones
2736 helper.SetSubShape( F );
2738 StdMeshers_FaceSide::GetFaceWires( F, *_mesh, /*skipMediumNodes=*/true,
2739 error, viscousMesh );
2740 error = AddSegmentsToMesh( *_ngMesh, occgeo, wires, helper, nodeVec );
2742 if ( !error ) error = SMESH_ComputeError::New();
2744 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2747 // Let netgen compute 2D mesh
2748 startWith = netgen::MESHCONST_MESHSURFACE;
2749 endWith = _optimize ? netgen::MESHCONST_OPTSURFACE : netgen::MESHCONST_MESHSURFACE;
2754 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2756 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2758 if(netgen::multithread.terminate)
2761 comment << text (err);
2763 catch (Standard_Failure& ex)
2765 comment << text(ex);
2766 //err = 1; -- try to make volumes anyway
2768 catch (netgen::NgException exc)
2770 comment << text(exc);
2771 //err = 1; -- try to make volumes anyway
2776 doneTime += faceMeshingTime + ( _optimize ? faceOptimizTime : 0 );
2777 _ticTime = doneTime / _totalTime / _progressTic;
2779 // ---------------------
2780 // generate volume mesh
2781 // ---------------------
2782 // Fill _ngMesh with nodes and faces of computed 2D submeshes
2783 if ( !err && _isVolume && ( !meshedSM[ MeshDim_2D ].empty() || mparams.quad ))
2785 // load SMESH with computed segments and faces
2786 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2788 // compute pyramids on quadrangles
2789 SMESH_ProxyMesh::Ptr proxyMesh;
2790 if ( _mesh->NbQuadrangles() > 0 )
2791 for ( int iS = 1; iS <= occgeo.somap.Extent(); ++iS )
2793 StdMeshers_QuadToTriaAdaptor* Adaptor = new StdMeshers_QuadToTriaAdaptor;
2794 proxyMesh.reset( Adaptor );
2796 int nbPyrams = _mesh->NbPyramids();
2797 Adaptor->Compute( *_mesh, occgeo.somap(iS) );
2798 if ( nbPyrams != _mesh->NbPyramids() )
2800 list< SMESH_subMesh* > quadFaceSM;
2801 for (TopExp_Explorer face(occgeo.somap(iS), TopAbs_FACE); face.More(); face.Next())
2802 if ( Adaptor->GetProxySubMesh( face.Current() ))
2804 quadFaceSM.push_back( _mesh->GetSubMesh( face.Current() ));
2805 meshedSM[ MeshDim_2D ].remove( quadFaceSM.back() );
2807 FillNgMesh(occgeo, *_ngMesh, nodeVec, quadFaceSM, &quadHelper, proxyMesh);
2810 // fill _ngMesh with faces of sub-meshes
2811 err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_2D ], &quadHelper));
2812 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2813 //toPython( _ngMesh, "/tmp/ngPython.py");
2815 if (!err && _isVolume)
2817 // Pass 3D simple parameters to NETGEN
2818 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
2819 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
2821 if ( double vol = simple3d->GetMaxElementVolume() ) {
2823 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
2824 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2827 // length from faces
2828 mparams.maxh = _ngMesh->AverageH();
2830 _ngMesh->SetGlobalH (mparams.maxh);
2831 mparams.grading = 0.4;
2833 _ngMesh->CalcLocalH(mparams.grading);
2835 _ngMesh->CalcLocalH();
2838 // Care of vertices internal in solids and internal faces (issue 0020676)
2839 if ( internals.hasInternalVertexInSolid() || internals.hasInternalFaces() )
2841 // store computed faces in SMESH in order not to create SMESH
2842 // faces for ng faces added here
2843 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2844 // add ng faces to solids with internal vertices
2845 AddIntVerticesInSolids( occgeo, *_ngMesh, nodeVec, internals );
2846 // duplicate mesh faces on internal faces
2847 FixIntFaces( occgeo, *_ngMesh, internals );
2848 initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2850 // Let netgen compute 3D mesh
2851 startWith = endWith = netgen::MESHCONST_MESHVOLUME;
2856 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2858 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2860 if(netgen::multithread.terminate)
2863 if ( comment.empty() ) // do not overwrite a previos error
2864 comment << text(err);
2866 catch (Standard_Failure& ex)
2868 if ( comment.empty() ) // do not overwrite a previos error
2869 comment << text(ex);
2872 catch (netgen::NgException exc)
2874 if ( comment.empty() ) // do not overwrite a previos error
2875 comment << text(exc);
2878 _ticTime = ( doneTime += voluMeshingTime ) / _totalTime / _progressTic;
2880 // Let netgen optimize 3D mesh
2881 if ( !err && _optimize )
2883 startWith = endWith = netgen::MESHCONST_OPTVOLUME;
2888 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2890 err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2892 if(netgen::multithread.terminate)
2895 if ( comment.empty() ) // do not overwrite a previos error
2896 comment << text(err);
2898 catch (Standard_Failure& ex)
2900 if ( comment.empty() ) // do not overwrite a previos error
2901 comment << text(ex);
2903 catch (netgen::NgException exc)
2905 if ( comment.empty() ) // do not overwrite a previos error
2906 comment << text(exc);
2910 if (!err && mparams.secondorder > 0)
2915 if ( !meshedSM[ MeshDim_1D ].empty() )
2917 // remove segments not attached to geometry (IPAL0052479)
2918 for (int i = 1; i <= _ngMesh->GetNSeg(); ++i)
2920 const netgen::Segment & seg = _ngMesh->LineSegment (i);
2921 if ( seg.epgeominfo[ 0 ].edgenr == 0 )
2922 _ngMesh->DeleteSegment( i );
2924 _ngMesh->Compress();
2926 // convert to quadratic
2927 netgen::OCCRefinementSurfaces ref (occgeo);
2928 ref.MakeSecondOrder (*_ngMesh);
2930 // care of elements already loaded to SMESH
2931 // if ( initState._nbSegments > 0 )
2932 // makeQuadratic( occgeo.emap, _mesh );
2933 // if ( initState._nbFaces > 0 )
2934 // makeQuadratic( occgeo.fmap, _mesh );
2936 catch (Standard_Failure& ex)
2938 if ( comment.empty() ) // do not overwrite a previos error
2939 comment << "Exception in netgen at passing to 2nd order ";
2941 catch (netgen::NgException exc)
2943 if ( comment.empty() ) // do not overwrite a previos error
2944 comment << exc.What();
2949 _ticTime = 0.98 / _progressTic;
2951 //int nbNod = _ngMesh->GetNP();
2952 //int nbSeg = _ngMesh->GetNSeg();
2953 int nbFac = _ngMesh->GetNSE();
2954 int nbVol = _ngMesh->GetNE();
2955 bool isOK = ( !err && (_isVolume ? (nbVol > 0) : (nbFac > 0)) );
2957 // Feed back the SMESHDS with the generated Nodes and Elements
2958 if ( true /*isOK*/ ) // get whatever built
2960 FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment, &quadHelper );
2962 if ( quadHelper.GetIsQuadratic() ) // remove free nodes
2963 for ( size_t i = 0; i < nodeVec.size(); ++i )
2964 if ( nodeVec[i] && nodeVec[i]->NbInverseElements() == 0 )
2965 _mesh->GetMeshDS()->RemoveFreeNode( nodeVec[i], 0, /*fromGroups=*/false );
2967 SMESH_ComputeErrorPtr readErr = ReadErrors(nodeVec);
2968 if ( readErr && !readErr->myBadElements.empty() )
2971 if ( !comment.empty() && !readErr->myComment.empty() ) comment += "\n";
2972 comment += readErr->myComment;
2974 if ( error->IsOK() && ( !isOK || comment.size() > 0 ))
2975 error->myName = COMPERR_ALGO_FAILED;
2976 if ( !comment.empty() )
2977 error->myComment = comment;
2979 // SetIsAlwaysComputed( true ) to empty sub-meshes, which
2980 // appear if the geometry contains coincident sub-shape due
2981 // to bool merge_solids = 1; in netgen/libsrc/occ/occgenmesh.cpp
2982 const int nbMaps = 2;
2983 const TopTools_IndexedMapOfShape* geoMaps[nbMaps] =
2984 { & occgeo.vmap, & occgeo.emap/*, & occgeo.fmap*/ };
2985 for ( int iMap = 0; iMap < nbMaps; ++iMap )
2986 for (int i = 1; i <= geoMaps[iMap]->Extent(); i++)
2987 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( geoMaps[iMap]->FindKey(i)))
2988 if ( !sm->IsMeshComputed() )
2989 sm->SetIsAlwaysComputed( true );
2991 // set bad compute error to subshapes of all failed sub-shapes
2992 if ( !error->IsOK() )
2994 bool pb2D = false, pb3D = false;
2995 for (int i = 1; i <= occgeo.fmap.Extent(); i++) {
2996 int status = occgeo.facemeshstatus[i-1];
2997 if (status == netgen::FACE_MESHED_OK ) continue;
2998 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.fmap( i ))) {
2999 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
3000 if ( !smError || smError->IsOK() ) {
3001 if ( status == netgen::FACE_FAILED )
3002 smError.reset( new SMESH_ComputeError( *error ));
3004 smError.reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED, "Ignored" ));
3005 if ( SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
3006 smError->myName = COMPERR_WARNING;
3008 pb2D = pb2D || smError->IsKO();
3011 if ( !pb2D ) // all faces are OK
3012 for (int i = 1; i <= occgeo.somap.Extent(); i++)
3013 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.somap( i )))
3015 bool smComputed = nbVol && !sm->IsEmpty();
3016 if ( smComputed && internals.hasInternalVertexInSolid( sm->GetId() ))
3018 int nbIntV = internals.getSolidsWithVertices().find( sm->GetId() )->second.size();
3019 SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
3020 smComputed = ( smDS->NbElements() > 0 || smDS->NbNodes() > nbIntV );
3022 SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
3023 if ( !smComputed && ( !smError || smError->IsOK() ))
3025 smError.reset( new SMESH_ComputeError( *error ));
3026 if ( nbVol && SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
3028 smError->myName = COMPERR_WARNING;
3030 else if ( !smError->myBadElements.empty() ) // bad surface mesh
3032 if ( !hasBadElemOnSolid( smError->myBadElements, sm ))
3036 pb3D = pb3D || ( smError && smError->IsKO() );
3038 if ( !pb2D && !pb3D )
3039 err = 0; // no fatal errors, only warnings
3042 ngLib._isComputeOk = !err;
3047 //=============================================================================
3051 //=============================================================================
3052 bool NETGENPlugin_Mesher::Evaluate(MapShapeNbElems& aResMap)
3054 netgen::MeshingParameters& mparams = netgen::mparam;
3057 // -------------------------
3058 // Prepare OCC geometry
3059 // -------------------------
3060 netgen::OCCGeometry occgeo;
3061 list< SMESH_subMesh* > meshedSM[4]; // for 0-3 dimensions
3062 NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
3063 PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
3065 bool tooManyElems = false;
3066 const int hugeNb = std::numeric_limits<int>::max() / 100;
3071 // pass 1D simple parameters to NETGEN
3074 // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
3075 mparams.uselocalh = false;
3076 mparams.grading = 0.8; // not limitited size growth
3078 if ( _simpleHyp->GetNumberOfSegments() )
3080 mparams.maxh = occgeo.boundingbox.Diam();
3083 mparams.maxh = _simpleHyp->GetLocalLength();
3086 if ( mparams.maxh == 0.0 )
3087 mparams.maxh = occgeo.boundingbox.Diam();
3088 if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
3089 mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
3091 // let netgen create _ngMesh and calculate element size on not meshed shapes
3092 NETGENPlugin_NetgenLibWrapper ngLib;
3093 netgen::Mesh *ngMesh = NULL;
3097 int startWith = netgen::MESHCONST_ANALYSE;
3098 int endWith = netgen::MESHCONST_MESHEDGES;
3100 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, mparams, startWith, endWith);
3102 int err = netgen::OCCGenerateMesh(occgeo, ngMesh, startWith, endWith, optstr);
3105 if(netgen::multithread.terminate)
3108 ngLib.setMesh(( Ng_Mesh*) ngMesh );
3110 if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( _shape ))
3111 sm->GetComputeError().reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED ));
3116 // Pass 1D simple parameters to NETGEN
3117 // --------------------------------
3118 int nbSeg = _simpleHyp->GetNumberOfSegments();
3119 double segSize = _simpleHyp->GetLocalLength();
3120 for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
3122 const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
3124 segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
3125 setLocalSize( e, segSize, *ngMesh );
3128 else // if ( ! _simpleHyp )
3130 // Local size on vertices and edges
3131 // --------------------------------
3132 for(std::map<int,double>::const_iterator it=EdgeId2LocalSize.begin(); it!=EdgeId2LocalSize.end(); it++)
3134 int key = (*it).first;
3135 double hi = (*it).second;
3136 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3137 const TopoDS_Edge& e = TopoDS::Edge(shape);
3138 setLocalSize( e, hi, *ngMesh );
3140 for(std::map<int,double>::const_iterator it=VertexId2LocalSize.begin(); it!=VertexId2LocalSize.end(); it++)
3142 int key = (*it).first;
3143 double hi = (*it).second;
3144 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3145 const TopoDS_Vertex& v = TopoDS::Vertex(shape);
3146 gp_Pnt p = BRep_Tool::Pnt(v);
3147 NETGENPlugin_Mesher::RestrictLocalSize( *ngMesh, p.XYZ(), hi );
3149 for(map<int,double>::const_iterator it=FaceId2LocalSize.begin();
3150 it!=FaceId2LocalSize.end(); it++)
3152 int key = (*it).first;
3153 double val = (*it).second;
3154 const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
3155 int faceNgID = occgeo.fmap.FindIndex(shape);
3156 occgeo.SetFaceMaxH(faceNgID, val);
3157 for ( TopExp_Explorer edgeExp( shape, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
3158 setLocalSize( TopoDS::Edge( edgeExp.Current() ), val, *ngMesh );
3161 // calculate total nb of segments and length of edges
3162 double fullLen = 0.0;
3164 int entity = mparams.secondorder > 0 ? SMDSEntity_Quad_Edge : SMDSEntity_Edge;
3165 TopTools_DataMapOfShapeInteger Edge2NbSeg;
3166 for (TopExp_Explorer exp(_shape, TopAbs_EDGE); exp.More(); exp.Next())
3168 TopoDS_Edge E = TopoDS::Edge( exp.Current() );
3169 if( !Edge2NbSeg.Bind(E,0) )
3172 double aLen = SMESH_Algo::EdgeLength(E);
3175 vector<int>& aVec = aResMap[_mesh->GetSubMesh(E)];
3177 aVec.resize( SMDSEntity_Last, 0);
3179 fullNbSeg += aVec[ entity ];
3182 // store nb of segments computed by Netgen
3183 NCollection_Map<Link> linkMap;
3184 for (int i = 1; i <= ngMesh->GetNSeg(); ++i )
3186 const netgen::Segment& seg = ngMesh->LineSegment(i);
3187 Link link(seg[0], seg[1]);
3188 if ( !linkMap.Add( link )) continue;
3189 int aGeomEdgeInd = seg.epgeominfo[0].edgenr;
3190 if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
3192 vector<int>& aVec = aResMap[_mesh->GetSubMesh(occgeo.emap(aGeomEdgeInd))];
3196 // store nb of nodes on edges computed by Netgen
3197 TopTools_DataMapIteratorOfDataMapOfShapeInteger Edge2NbSegIt(Edge2NbSeg);
3198 for (; Edge2NbSegIt.More(); Edge2NbSegIt.Next())
3200 vector<int>& aVec = aResMap[_mesh->GetSubMesh(Edge2NbSegIt.Key())];
3201 if ( aVec[ entity ] > 1 && aVec[ SMDSEntity_Node ] == 0 )
3202 aVec[SMDSEntity_Node] = mparams.secondorder > 0 ? 2*aVec[ entity ]-1 : aVec[ entity ]-1;
3204 fullNbSeg += aVec[ entity ];
3205 Edge2NbSeg( Edge2NbSegIt.Key() ) = aVec[ entity ];
3207 if ( fullNbSeg == 0 )
3214 if ( double area = _simpleHyp->GetMaxElementArea() ) {
3216 mparams.maxh = sqrt(2. * area/sqrt(3.0));
3217 mparams.grading = 0.4; // moderate size growth
3220 // length from edges
3221 mparams.maxh = fullLen/fullNbSeg;
3222 mparams.grading = 0.2; // slow size growth
3225 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3226 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3228 for (TopExp_Explorer exp(_shape, TopAbs_FACE); exp.More(); exp.Next())
3230 TopoDS_Face F = TopoDS::Face( exp.Current() );
3231 SMESH_subMesh *sm = _mesh->GetSubMesh(F);
3233 BRepGProp::SurfaceProperties(F,G);
3234 double anArea = G.Mass();
3235 tooManyElems = tooManyElems || ( anArea/hugeNb > mparams.maxh*mparams.maxh );
3237 if ( !tooManyElems )
3239 TopTools_MapOfShape egdes;
3240 for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next())
3241 if ( egdes.Add( exp1.Current() ))
3242 nb1d += Edge2NbSeg.Find(exp1.Current());
3244 int nbFaces = tooManyElems ? hugeNb : int( 4*anArea / (mparams.maxh*mparams.maxh*sqrt(3.)));
3245 int nbNodes = tooManyElems ? hugeNb : (( nbFaces*3 - (nb1d-1)*2 ) / 6 + 1 );
3247 vector<int> aVec(SMDSEntity_Last, 0);
3248 if( mparams.secondorder > 0 ) {
3249 int nb1d_in = (nbFaces*3 - nb1d) / 2;
3250 aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
3251 aVec[SMDSEntity_Quad_Triangle] = nbFaces;
3254 aVec[SMDSEntity_Node] = Max ( nbNodes, 0 );
3255 aVec[SMDSEntity_Triangle] = nbFaces;
3257 aResMap[sm].swap(aVec);
3264 // pass 3D simple parameters to NETGEN
3265 const NETGENPlugin_SimpleHypothesis_3D* simple3d =
3266 dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
3268 if ( double vol = simple3d->GetMaxElementVolume() ) {
3270 mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
3271 mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
3274 // using previous length from faces
3276 mparams.grading = 0.4;
3277 mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
3280 BRepGProp::VolumeProperties(_shape,G);
3281 double aVolume = G.Mass();
3282 double tetrVol = 0.1179*mparams.maxh*mparams.maxh*mparams.maxh;
3283 tooManyElems = tooManyElems || ( aVolume/hugeNb > tetrVol );
3284 int nbVols = tooManyElems ? hugeNb : int(aVolume/tetrVol);
3285 int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
3286 vector<int> aVec(SMDSEntity_Last, 0 );
3287 if ( tooManyElems ) // avoid FPE
3289 aVec[SMDSEntity_Node] = hugeNb;
3290 aVec[ mparams.secondorder > 0 ? SMDSEntity_Quad_Tetra : SMDSEntity_Tetra] = hugeNb;
3294 if( mparams.secondorder > 0 ) {
3295 aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
3296 aVec[SMDSEntity_Quad_Tetra] = nbVols;
3299 aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
3300 aVec[SMDSEntity_Tetra] = nbVols;
3303 SMESH_subMesh *sm = _mesh->GetSubMesh(_shape);
3304 aResMap[sm].swap(aVec);
3310 double NETGENPlugin_Mesher::GetProgress(const SMESH_Algo* holder,
3311 const int * algoProgressTic,
3312 const double * algoProgress) const
3314 ((int&) _progressTic ) = *algoProgressTic + 1;
3316 if ( !_occgeom ) return 0;
3318 double progress = -1;
3321 if ( _ticTime < 0 && netgen::multithread.task[0] == 'O'/*Optimizing surface*/ )
3323 ((double&) _ticTime ) = edgeFaceMeshingTime / _totalTime / _progressTic;
3325 else if ( !_optimize /*&& _occgeom->fmap.Extent() > 1*/ )
3327 int doneShapeIndex = -1;
3328 while ( doneShapeIndex+1 < _occgeom->facemeshstatus.Size() &&
3329 _occgeom->facemeshstatus[ doneShapeIndex+1 ])
3331 if ( doneShapeIndex+1 != _curShapeIndex )
3333 ((int&) _curShapeIndex) = doneShapeIndex+1;
3334 double doneShapeRate = _curShapeIndex / double( _occgeom->fmap.Extent() );
3335 double doneTime = edgeMeshingTime + doneShapeRate * faceMeshingTime;
3336 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3337 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3338 // << " " << doneTime / _totalTime / _progressTic << endl;
3342 else if ( !_optimize && _occgeom->somap.Extent() > 1 )
3344 int curShapeIndex = _curShapeIndex;
3345 if ( _ngMesh->GetNE() > 0 )
3347 netgen::Element el = (*_ngMesh)[netgen::ElementIndex( _ngMesh->GetNE()-1 )];
3348 curShapeIndex = el.GetIndex();
3350 if ( curShapeIndex != _curShapeIndex )
3352 ((int&) _curShapeIndex) = curShapeIndex;
3353 double doneShapeRate = _curShapeIndex / double( _occgeom->somap.Extent() );
3354 double doneTime = edgeFaceMeshingTime + doneShapeRate * voluMeshingTime;
3355 ((double&) _ticTime) = doneTime / _totalTime / _progressTic;
3356 // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3357 // << " " << doneTime / _totalTime / _progressTic << endl;
3361 progress = Max( *algoProgressTic * _ticTime, *algoProgress );
3364 ((int&) *algoProgressTic )++;
3365 ((double&) *algoProgress) = progress;
3367 //cout << progress << " " << *algoProgressTic << " " << netgen::multithread.task << " "<< _ticTime << endl;
3369 return Min( progress, 0.99 );
3372 //================================================================================
3374 * \brief Read mesh entities preventing successful computation from "test.out" file
3376 //================================================================================
3378 SMESH_ComputeErrorPtr
3379 NETGENPlugin_Mesher::ReadErrors(const vector<const SMDS_MeshNode* >& nodeVec)
3381 SMESH_ComputeErrorPtr err = SMESH_ComputeError::New
3382 (COMPERR_BAD_INPUT_MESH, "Some edges multiple times in surface mesh");
3383 SMESH_File file("test.out");
3385 vector<int> three1(3), three2(3);
3386 const char* badEdgeStr = " multiple times in surface mesh";
3387 const int badEdgeStrLen = strlen( badEdgeStr );
3388 const int nbNodes = nodeVec.size();
3390 while( !file.eof() )
3392 if ( strncmp( file, "Edge ", 5 ) == 0 &&
3393 file.getInts( two ) &&
3394 strncmp( file, badEdgeStr, badEdgeStrLen ) == 0 &&
3395 two[0] < nbNodes && two[1] < nbNodes )
3397 err->myBadElements.push_back( new SMDS_LinearEdge( nodeVec[ two[0]], nodeVec[ two[1]] ));
3398 file += badEdgeStrLen;
3400 else if ( strncmp( file, "Intersecting: ", 14 ) == 0 )
3403 // openelement 18 with open element 126
3407 const char* pos = file;
3408 bool ok = ( strncmp( file, "openelement ", 12 ) == 0 );
3409 ok = ok && file.getInts( two );
3410 ok = ok && file.getInts( three1 );
3411 ok = ok && file.getInts( three2 );
3412 for ( int i = 0; ok && i < 3; ++i )
3413 ok = ( three1[i] < nbNodes && nodeVec[ three1[i]]);
3414 for ( int i = 0; ok && i < 3; ++i )
3415 ok = ( three2[i] < nbNodes && nodeVec[ three2[i]]);
3418 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three1[0]],
3419 nodeVec[ three1[1]],
3420 nodeVec[ three1[2]]));
3421 err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three2[0]],
3422 nodeVec[ three2[1]],
3423 nodeVec[ three2[2]]));
3424 err->myComment = "Intersecting triangles";
3438 size_t nbBadElems = err->myBadElements.size();
3439 if ( nbBadElems ) nbBadElems++; // avoid warning: variable set but not used
3445 //================================================================================
3447 * \brief Write a python script creating an equivalent SALOME mesh.
3448 * This is useful to see what mesh is passed as input for the next step of mesh
3449 * generation (of mesh of higher dimension)
3451 //================================================================================
3453 void NETGENPlugin_Mesher::toPython( const netgen::Mesh* ngMesh,
3454 const std::string& pyFile)
3456 ofstream outfile(pyFile.c_str(), ios::out);
3457 if ( !outfile ) return;
3459 outfile << "import SMESH" << endl
3460 << "from salome.smesh import smeshBuilder" << endl
3461 << "smesh = smeshBuilder.New(salome.myStudy)" << endl
3462 << "mesh = smesh.Mesh()" << endl << endl;
3464 using namespace netgen;
3466 for (pi = PointIndex::BASE;
3467 pi < ngMesh->GetNP()+PointIndex::BASE; pi++)
3469 outfile << "mesh.AddNode( ";
3470 outfile << (*ngMesh)[pi](0) << ", ";
3471 outfile << (*ngMesh)[pi](1) << ", ";
3472 outfile << (*ngMesh)[pi](2) << ") ## "<< pi << endl;
3475 int nbDom = ngMesh->GetNDomains();
3476 for ( int i = 0; i < nbDom; ++i )
3477 outfile<< "grp" << i+1 << " = mesh.CreateEmptyGroup( SMESH.FACE, 'domain"<< i+1 << "')"<< endl;
3479 SurfaceElementIndex sei;
3480 for (sei = 0; sei < ngMesh->GetNSE(); sei++)
3482 outfile << "mesh.AddFace([ ";
3483 Element2d sel = (*ngMesh)[sei];
3484 for (int j = 0; j < sel.GetNP(); j++)
3485 outfile << sel[j] << ( j+1 < sel.GetNP() ? ", " : " ])");
3486 if ( sel.IsDeleted() ) outfile << " ## IsDeleted ";
3489 if ((*ngMesh)[sei].GetIndex())
3491 if ( int dom1 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainIn())
3492 outfile << "grp"<< dom1 <<".Add([ " << (int)sei+1 << " ])" << endl;
3493 if ( int dom2 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainOut())
3494 outfile << "grp"<< dom2 <<".Add([ " << (int)sei+1 << " ])" << endl;
3498 for (ElementIndex ei = 0; ei < ngMesh->GetNE(); ei++)
3500 Element el = (*ngMesh)[ei];
3501 outfile << "mesh.AddVolume([ ";
3502 for (int j = 0; j < el.GetNP(); j++)
3503 outfile << el[j] << ( j+1 < el.GetNP() ? ", " : " ])");
3507 for (int i = 1; i <= ngMesh->GetNSeg(); i++)
3509 const Segment & seg = ngMesh->LineSegment (i);
3510 outfile << "mesh.AddEdge([ "
3512 << seg[1] << " ])" << endl;
3514 cout << "Write " << pyFile << endl;
3517 //================================================================================
3519 * \brief Constructor of NETGENPlugin_ngMeshInfo
3521 //================================================================================
3523 NETGENPlugin_ngMeshInfo::NETGENPlugin_ngMeshInfo( netgen::Mesh* ngMesh):
3528 _nbNodes = ngMesh->GetNP();
3529 _nbSegments = ngMesh->GetNSeg();
3530 _nbFaces = ngMesh->GetNSE();
3531 _nbVolumes = ngMesh->GetNE();
3535 _nbNodes = _nbSegments = _nbFaces = _nbVolumes = 0;
3539 //================================================================================
3541 * \brief Copy LocalH member from one netgen mesh to another
3543 //================================================================================
3545 void NETGENPlugin_ngMeshInfo::transferLocalH( netgen::Mesh* fromMesh,
3546 netgen::Mesh* toMesh )
3548 if ( !fromMesh->LocalHFunctionGenerated() ) return;
3549 if ( !toMesh->LocalHFunctionGenerated() )
3551 toMesh->CalcLocalH(netgen::mparam.grading);
3553 toMesh->CalcLocalH();
3556 const size_t size = sizeof( netgen::LocalH );
3557 _copyOfLocalH = new char[ size ];
3558 memcpy( (void*)_copyOfLocalH, (void*)&toMesh->LocalHFunction(), size );
3559 memcpy( (void*)&toMesh->LocalHFunction(), (void*)&fromMesh->LocalHFunction(), size );
3562 //================================================================================
3564 * \brief Restore LocalH member of a netgen mesh
3566 //================================================================================
3568 void NETGENPlugin_ngMeshInfo::restoreLocalH( netgen::Mesh* toMesh )
3570 if ( _copyOfLocalH )
3572 const size_t size = sizeof( netgen::LocalH );
3573 memcpy( (void*)&toMesh->LocalHFunction(), (void*)_copyOfLocalH, size );
3574 delete [] _copyOfLocalH;
3579 //================================================================================
3581 * \brief Find "internal" sub-shapes
3583 //================================================================================
3585 NETGENPlugin_Internals::NETGENPlugin_Internals( SMESH_Mesh& mesh,
3586 const TopoDS_Shape& shape,
3588 : _mesh( mesh ), _is3D( is3D )
3590 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3592 TopExp_Explorer f,e;
3593 for ( f.Init( shape, TopAbs_FACE ); f.More(); f.Next() )
3595 int faceID = meshDS->ShapeToIndex( f.Current() );
3597 // find not computed internal edges
3599 for ( e.Init( f.Current().Oriented(TopAbs_FORWARD), TopAbs_EDGE ); e.More(); e.Next() )
3600 if ( e.Current().Orientation() == TopAbs_INTERNAL )
3602 SMESH_subMesh* eSM = mesh.GetSubMesh( e.Current() );
3603 if ( eSM->IsEmpty() )
3605 _e2face.insert( make_pair( eSM->GetId(), faceID ));
3606 for ( TopoDS_Iterator v(e.Current()); v.More(); v.Next() )
3607 _e2face.insert( make_pair( meshDS->ShapeToIndex( v.Value() ), faceID ));
3611 // find internal vertices in a face
3612 set<int> intVV; // issue 0020850 where same vertex is twice in a face
3613 for ( TopoDS_Iterator fSub( f.Current() ); fSub.More(); fSub.Next())
3614 if ( fSub.Value().ShapeType() == TopAbs_VERTEX )
3616 int vID = meshDS->ShapeToIndex( fSub.Value() );
3617 if ( intVV.insert( vID ).second )
3618 _f2v[ faceID ].push_back( vID );
3623 // find internal faces and their subshapes where nodes are to be doubled
3624 // to make a crack with non-sewed borders
3626 if ( f.Current().Orientation() == TopAbs_INTERNAL )
3628 _intShapes.insert( meshDS->ShapeToIndex( f.Current() ));
3631 list< TopoDS_Shape > edges;
3632 for ( e.Init( f.Current(), TopAbs_EDGE ); e.More(); e.Next())
3633 if ( SMESH_MesherHelper::NbAncestors( e.Current(), mesh, TopAbs_FACE ) > 1 )
3635 _intShapes.insert( meshDS->ShapeToIndex( e.Current() ));
3636 edges.push_back( e.Current() );
3637 // find border faces
3638 PShapeIteratorPtr fIt =
3639 SMESH_MesherHelper::GetAncestors( edges.back(),mesh,TopAbs_FACE );
3640 while ( const TopoDS_Shape* pFace = fIt->next() )
3641 if ( !pFace->IsSame( f.Current() ))
3642 _borderFaces.insert( meshDS->ShapeToIndex( *pFace ));
3645 // we consider vertex internal if it is shared by more than one internal edge
3646 list< TopoDS_Shape >::iterator edge = edges.begin();
3647 for ( ; edge != edges.end(); ++edge )
3648 for ( TopoDS_Iterator v( *edge ); v.More(); v.Next() )
3650 set<int> internalEdges;
3651 PShapeIteratorPtr eIt =
3652 SMESH_MesherHelper::GetAncestors( v.Value(),mesh,TopAbs_EDGE );
3653 while ( const TopoDS_Shape* pEdge = eIt->next() )
3655 int edgeID = meshDS->ShapeToIndex( *pEdge );
3656 if ( isInternalShape( edgeID ))
3657 internalEdges.insert( edgeID );
3659 if ( internalEdges.size() > 1 )
3660 _intShapes.insert( meshDS->ShapeToIndex( v.Value() ));
3664 } // loop on geom faces
3666 // find vertices internal in solids
3669 for ( TopExp_Explorer so(shape, TopAbs_SOLID); so.More(); so.Next())
3671 int soID = meshDS->ShapeToIndex( so.Current() );
3672 for ( TopoDS_Iterator soSub( so.Current() ); soSub.More(); soSub.Next())
3673 if ( soSub.Value().ShapeType() == TopAbs_VERTEX )
3674 _s2v[ soID ].push_back( meshDS->ShapeToIndex( soSub.Value() ));
3679 //================================================================================
3681 * \brief Find mesh faces on non-internal geom faces sharing internal edge
3682 * some nodes of which are to be doubled to make the second border of the "crack"
3684 //================================================================================
3686 void NETGENPlugin_Internals::findBorderElements( TIDSortedElemSet & borderElems )
3688 if ( _intShapes.empty() ) return;
3690 SMESH_Mesh& mesh = const_cast<SMESH_Mesh&>(_mesh);
3691 SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3693 // loop on internal geom edges
3694 set<int>::const_iterator intShapeId = _intShapes.begin();
3695 for ( ; intShapeId != _intShapes.end(); ++intShapeId )
3697 const TopoDS_Shape& s = meshDS->IndexToShape( *intShapeId );
3698 if ( s.ShapeType() != TopAbs_EDGE ) continue;
3700 // get internal and non-internal geom faces sharing the internal edge <s>
3702 set<int>::iterator bordFace = _borderFaces.end();
3703 PShapeIteratorPtr faces = SMESH_MesherHelper::GetAncestors( s, _mesh, TopAbs_FACE );
3704 while ( const TopoDS_Shape* pFace = faces->next() )
3706 int faceID = meshDS->ShapeToIndex( *pFace );
3707 if ( isInternalShape( faceID ))
3710 bordFace = _borderFaces.insert( faceID ).first;
3712 if ( bordFace == _borderFaces.end() || !intFace ) continue;
3714 // get all links of mesh faces on internal geom face sharing nodes on edge <s>
3715 set< SMESH_OrientedLink > links; //!< links of faces on internal geom face
3716 list<const SMDS_MeshElement*> suspectFaces[2]; //!< mesh faces on border geom faces
3717 int nbSuspectFaces = 0;
3718 SMESHDS_SubMesh* intFaceSM = meshDS->MeshElements( intFace );
3719 if ( !intFaceSM || intFaceSM->NbElements() == 0 ) continue;
3720 SMESH_subMeshIteratorPtr smIt = mesh.GetSubMesh( s )->getDependsOnIterator(true,true);
3721 while ( smIt->more() )
3723 SMESHDS_SubMesh* sm = smIt->next()->GetSubMeshDS();
3724 if ( !sm ) continue;
3725 SMDS_NodeIteratorPtr nIt = sm->GetNodes();
3726 while ( nIt->more() )
3728 const SMDS_MeshNode* nOnEdge = nIt->next();
3729 SMDS_ElemIteratorPtr fIt = nOnEdge->GetInverseElementIterator(SMDSAbs_Face);
3730 while ( fIt->more() )
3732 const SMDS_MeshElement* f = fIt->next();
3733 const int nbNodes = f->NbCornerNodes();
3734 if ( intFaceSM->Contains( f ))
3736 for ( int i = 0; i < nbNodes; ++i )
3737 links.insert( SMESH_OrientedLink( f->GetNode(i), f->GetNode((i+1)%nbNodes)));
3742 for ( int i = 0; i < nbNodes; ++i )
3743 nbDblNodes += isInternalShape( f->GetNode(i)->getshapeId() );
3745 suspectFaces[ nbDblNodes < 2 ].push_back( f );
3751 // suspectFaces[0] having link with same orientation as mesh faces on
3752 // the internal geom face are <borderElems>. suspectFaces[1] have
3753 // only one node on edge <s>, we decide on them later (at the 2nd loop)
3754 // by links of <borderElems> found at the 1st and 2nd loops
3755 set< SMESH_OrientedLink > borderLinks;
3756 for ( int isPostponed = 0; isPostponed < 2; ++isPostponed )
3758 list<const SMDS_MeshElement*>::iterator fIt = suspectFaces[isPostponed].begin();
3759 for ( int nbF = 0; fIt != suspectFaces[isPostponed].end(); ++fIt, ++nbF )
3761 const SMDS_MeshElement* f = *fIt;
3762 bool isBorder = false, linkFound = false, borderLinkFound = false;
3763 list< SMESH_OrientedLink > faceLinks;
3764 int nbNodes = f->NbCornerNodes();
3765 for ( int i = 0; i < nbNodes; ++i )
3767 SMESH_OrientedLink link( f->GetNode(i), f->GetNode((i+1)%nbNodes));
3768 faceLinks.push_back( link );
3771 set< SMESH_OrientedLink >::iterator foundLink = links.find( link );
3772 if ( foundLink != links.end() )
3775 isBorder = ( foundLink->_reversed == link._reversed );
3776 if ( !isBorder && !isPostponed ) break;
3777 faceLinks.pop_back();
3779 else if ( isPostponed && !borderLinkFound )
3781 foundLink = borderLinks.find( link );
3782 if ( foundLink != borderLinks.end() )
3784 borderLinkFound = true;
3785 isBorder = ( foundLink->_reversed != link._reversed );
3792 borderElems.insert( f );
3793 borderLinks.insert( faceLinks.begin(), faceLinks.end() );
3795 else if ( !linkFound && !borderLinkFound )
3797 suspectFaces[1].push_back( f );
3798 if ( nbF > 2 * nbSuspectFaces )
3799 break; // dead loop protection
3806 //================================================================================
3808 * \brief put internal shapes in maps and fill in submeshes to precompute
3810 //================================================================================
3812 void NETGENPlugin_Internals::getInternalEdges( TopTools_IndexedMapOfShape& fmap,
3813 TopTools_IndexedMapOfShape& emap,
3814 TopTools_IndexedMapOfShape& vmap,
3815 list< SMESH_subMesh* > smToPrecompute[])
3817 if ( !hasInternalEdges() ) return;
3818 map<int,int>::const_iterator ev_face = _e2face.begin();
3819 for ( ; ev_face != _e2face.end(); ++ev_face )
3821 const TopoDS_Shape& ev = _mesh.GetMeshDS()->IndexToShape( ev_face->first );
3822 const TopoDS_Shape& face = _mesh.GetMeshDS()->IndexToShape( ev_face->second );
3824 ( ev.ShapeType() == TopAbs_EDGE ? emap : vmap ).Add( ev );
3826 //cout<<"INTERNAL EDGE or VERTEX "<<ev_face->first<<" on face "<<ev_face->second<<endl;
3828 smToPrecompute[ MeshDim_1D ].push_back( _mesh.GetSubMeshContaining( ev_face->first ));
3832 //================================================================================
3834 * \brief return shapes and submeshes to be meshed and already meshed boundary submeshes
3836 //================================================================================
3838 void NETGENPlugin_Internals::getInternalFaces( TopTools_IndexedMapOfShape& fmap,
3839 TopTools_IndexedMapOfShape& emap,
3840 list< SMESH_subMesh* >& intFaceSM,
3841 list< SMESH_subMesh* >& boundarySM)
3843 if ( !hasInternalFaces() ) return;
3845 // <fmap> and <emap> are for not yet meshed shapes
3846 // <intFaceSM> is for submeshes of faces
3847 // <boundarySM> is for meshed edges and vertices
3852 set<int> shapeIDs ( _intShapes );
3853 if ( !_borderFaces.empty() )
3854 shapeIDs.insert( _borderFaces.begin(), _borderFaces.end() );
3856 set<int>::const_iterator intS = shapeIDs.begin();
3857 for ( ; intS != shapeIDs.end(); ++intS )
3859 SMESH_subMesh* sm = _mesh.GetSubMeshContaining( *intS );
3861 if ( sm->GetSubShape().ShapeType() != TopAbs_FACE ) continue;
3863 intFaceSM.push_back( sm );
3865 // add submeshes of not computed internal faces
3866 if ( !sm->IsEmpty() ) continue;
3868 SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(true,true);
3869 while ( smIt->more() )
3872 const TopoDS_Shape& s = sm->GetSubShape();
3874 if ( sm->IsEmpty() )
3877 switch ( s.ShapeType() ) {
3878 case TopAbs_FACE: fmap.Add ( s ); break;
3879 case TopAbs_EDGE: emap.Add ( s ); break;
3885 if ( s.ShapeType() != TopAbs_FACE )
3886 boundarySM.push_back( sm );
3892 //================================================================================
3894 * \brief Return true if given shape is to be precomputed in order to be correctly
3895 * added to netgen mesh
3897 //================================================================================
3899 bool NETGENPlugin_Internals::isShapeToPrecompute(const TopoDS_Shape& s)
3901 int shapeID = _mesh.GetMeshDS()->ShapeToIndex( s );
3902 switch ( s.ShapeType() ) {
3903 case TopAbs_FACE : break; //return isInternalShape( shapeID ) || isBorderFace( shapeID );
3904 case TopAbs_EDGE : return isInternalEdge( shapeID );
3905 case TopAbs_VERTEX: break;
3911 //================================================================================
3913 * \brief Return SMESH
3915 //================================================================================
3917 SMESH_Mesh& NETGENPlugin_Internals::getMesh() const
3919 return const_cast<SMESH_Mesh&>( _mesh );
3922 //================================================================================
3924 * \brief Access to a counter of NETGENPlugin_NetgenLibWrapper instances
3926 //================================================================================
3928 int& NETGENPlugin_NetgenLibWrapper::instanceCounter()
3930 static int theCouner = 0;
3934 //================================================================================
3936 * \brief Initialize netgen library
3938 //================================================================================
3940 NETGENPlugin_NetgenLibWrapper::NETGENPlugin_NetgenLibWrapper()
3942 if ( instanceCounter() == 0 )
3945 ++instanceCounter();
3947 _isComputeOk = false;
3951 if ( !getenv( "KEEP_NETGEN_OUTPUT" ))
3953 // redirect all netgen output (mycout,myerr,cout) to _outputFileName
3954 _outputFileName = getOutputFileName();
3955 _ngcout = netgen::mycout;
3956 _ngcerr = netgen::myerr;
3957 netgen::mycout = new ofstream ( _outputFileName.c_str() );
3958 netgen::myerr = netgen::mycout;
3959 _coutBuffer = std::cout.rdbuf();
3961 cout << "NOTE: netgen output is redirected to file " << _outputFileName << endl;
3963 std::cout.rdbuf( netgen::mycout->rdbuf() );
3967 _ngMesh = Ng_NewMesh();
3970 //================================================================================
3972 * \brief Finish using netgen library
3974 //================================================================================
3976 NETGENPlugin_NetgenLibWrapper::~NETGENPlugin_NetgenLibWrapper()
3978 --instanceCounter();
3980 Ng_DeleteMesh( _ngMesh );
3984 std::cout.rdbuf( _coutBuffer );
3991 //================================================================================
3993 * \brief Set netgen mesh to delete at destruction
3995 //================================================================================
3997 void NETGENPlugin_NetgenLibWrapper::setMesh( Ng_Mesh* mesh )
4000 Ng_DeleteMesh( _ngMesh );
4004 //================================================================================
4006 * \brief Return a unique file name
4008 //================================================================================
4010 std::string NETGENPlugin_NetgenLibWrapper::getOutputFileName()
4012 std::string aTmpDir = SALOMEDS_Tool::GetTmpDir();
4014 TCollection_AsciiString aGenericName = (char*)aTmpDir.c_str();
4015 aGenericName += "NETGEN_";
4017 aGenericName += getpid();
4019 aGenericName += _getpid();
4021 aGenericName += "_";
4022 aGenericName += Abs((Standard_Integer)(long) aGenericName.ToCString());
4023 aGenericName += ".out";
4025 return aGenericName.ToCString();
4028 //================================================================================
4030 * \brief Remove "test.out" and "problemfaces" files in current directory
4032 //================================================================================
4034 void NETGENPlugin_NetgenLibWrapper::RemoveTmpFiles()
4036 bool rm = SMESH_File("test.out").remove() ;
4038 if ( rm && netgen::testout && instanceCounter() == 0 )
4040 delete netgen::testout;
4041 netgen::testout = 0;
4044 SMESH_File("problemfaces").remove();
4045 SMESH_File("occmesh.rep").remove();
4048 //================================================================================
4050 * \brief Remove file with netgen output
4052 //================================================================================
4054 void NETGENPlugin_NetgenLibWrapper::removeOutputFile()
4056 if ( !_outputFileName.empty() )
4060 delete netgen::mycout;
4061 netgen::mycout = _ngcout;
4062 netgen::myerr = _ngcerr;
4065 string tmpDir = SALOMEDS_Tool::GetDirFromPath ( _outputFileName );
4066 string aFileName = SALOMEDS_Tool::GetNameFromPath( _outputFileName ) + ".out";
4067 SALOMEDS::ListOfFileNames_var aFiles = new SALOMEDS::ListOfFileNames;
4069 aFiles[0] = aFileName.c_str();
4071 SALOMEDS_Tool::RemoveTemporaryFiles( tmpDir.c_str(), aFiles.in(), true );