]> SALOME platform Git repositories - plugins/netgenplugin.git/blob - src/NETGENPlugin/NETGENPlugin_Mesher.cxx
Salome HOME
-#if (OCC_VERSION_MAJOR << 16 | OCC_VERSION_MINOR << 8 | OCC_VERSION_MAINTENANCE...
[plugins/netgenplugin.git] / src / NETGENPlugin / NETGENPlugin_Mesher.cxx
1 // Copyright (C) 2007-2013  CEA/DEN, EDF R&D, OPEN CASCADE
2 //
3 // Copyright (C) 2003-2007  OPEN CASCADE, EADS/CCR, LIP6, CEA/DEN,
4 // CEDRAT, EDF R&D, LEG, PRINCIPIA R&D, BUREAU VERITAS
5 //
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.
10 //
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.
15 //
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
19 //
20 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
21 //
22
23 //  NETGENPlugin : C++ implementation
24 // File      : NETGENPlugin_Mesher.cxx
25 // Author    : Michael Sazonov (OCN)
26 // Date      : 31/03/2006
27 // Project   : SALOME
28 //=============================================================================
29
30 #include "NETGENPlugin_Mesher.hxx"
31 #include "NETGENPlugin_Hypothesis_2D.hxx"
32 #include "NETGENPlugin_SimpleHypothesis_3D.hxx"
33
34 #include <SMDS_FaceOfNodes.hxx>
35 #include <SMDS_MeshElement.hxx>
36 #include <SMDS_MeshNode.hxx>
37 #include <SMESHDS_Mesh.hxx>
38 #include <SMESH_Block.hxx>
39 #include <SMESH_Comment.hxx>
40 #include <SMESH_ComputeError.hxx>
41 #include <SMESH_File.hxx>
42 #include <SMESH_Gen_i.hxx>
43 #include <SMESH_Mesh.hxx>
44 #include <SMESH_MesherHelper.hxx>
45 #include <SMESH_subMesh.hxx>
46 #include <StdMeshers_QuadToTriaAdaptor.hxx>
47 #include <StdMeshers_ViscousLayers2D.hxx>
48
49 #include <SALOMEDS_Tool.hxx>
50
51 #include <utilities.h>
52
53 #include <BRepBuilderAPI_Copy.hxx>
54 #include <BRep_Tool.hxx>
55 #include <Bnd_B3d.hxx>
56 #include <NCollection_Map.hxx>
57 #include <Standard_ErrorHandler.hxx>
58 #include <Standard_ProgramError.hxx>
59 #include <TopExp.hxx>
60 #include <TopExp_Explorer.hxx>
61 #include <TopTools_DataMapIteratorOfDataMapOfShapeInteger.hxx>
62 #include <TopTools_DataMapIteratorOfDataMapOfShapeShape.hxx>
63 #include <TopTools_DataMapOfShapeInteger.hxx>
64 #include <TopTools_DataMapOfShapeShape.hxx>
65 #include <TopTools_MapOfShape.hxx>
66 #include <TopoDS.hxx>
67 #include <OSD_File.hxx>
68 #include <OSD_Path.hxx>
69
70 // Netgen include files
71 #ifndef OCCGEOMETRY
72 #define OCCGEOMETRY
73 #endif
74 #include <occgeom.hpp>
75 #include <meshing.hpp>
76 //#include <ngexception.hpp>
77 namespace netgen {
78 #ifdef NETGEN_V5
79   extern int OCCGenerateMesh (OCCGeometry&, Mesh*&, MeshingParameters&, int, int);
80 #else
81   extern int OCCGenerateMesh (OCCGeometry&, Mesh*&, int, int, char*);
82 #endif
83   //extern void OCCSetLocalMeshSize(OCCGeometry & geom, Mesh & mesh);
84   extern MeshingParameters mparam;
85   extern volatile multithreadt multithread;
86 }
87
88 #include <vector>
89 #include <limits>
90
91 #ifdef WIN32
92 #include <process.h>
93 #endif
94 using namespace nglib;
95 using namespace std;
96
97 #ifdef _DEBUG_
98 #define nodeVec_ACCESS(index) ((SMDS_MeshNode*) nodeVec.at((index)))
99 #else
100 #define nodeVec_ACCESS(index) ((SMDS_MeshNode*) nodeVec[index])
101 #endif
102
103 #define NGPOINT_COORDS(p) p(0),p(1),p(2)
104
105 // dump elements added to ng mesh
106 //#define DUMP_SEGMENTS
107 //#define DUMP_TRIANGLES
108 //#define DUMP_TRIANGLES_SCRIPT "/tmp/trias.py" //!< debug AddIntVerticesInSolids()
109
110 TopTools_IndexedMapOfShape ShapesWithLocalSize;
111 std::map<int,double> VertexId2LocalSize;
112 std::map<int,double> EdgeId2LocalSize;
113 std::map<int,double> FaceId2LocalSize;
114
115 //=============================================================================
116 /*!
117  *
118  */
119 //=============================================================================
120
121 NETGENPlugin_Mesher::NETGENPlugin_Mesher (SMESH_Mesh*         mesh,
122                                           const TopoDS_Shape& aShape,
123                                           const bool          isVolume)
124   : _mesh    (mesh),
125     _shape   (aShape),
126     _isVolume(isVolume),
127     _optimize(true),
128     _fineness(NETGENPlugin_Hypothesis::GetDefaultFineness()),
129     _isViscousLayers2D(false),
130     _ngMesh(NULL),
131     _occgeom(NULL),
132     _curShapeIndex(-1),
133     _progressTic(1),
134     _totalTime(1.0),
135     _simpleHyp(NULL),
136     _ptrToMe(NULL)
137 {
138   SetDefaultParameters();
139   ShapesWithLocalSize.Clear();
140   VertexId2LocalSize.clear();
141   EdgeId2LocalSize.clear();
142   FaceId2LocalSize.clear();
143 }
144
145 //================================================================================
146 /*!
147  * Destuctor
148  */
149 //================================================================================
150
151 NETGENPlugin_Mesher::~NETGENPlugin_Mesher()
152 {
153   if ( _ptrToMe )
154     *_ptrToMe = NULL;
155   _ptrToMe = 0;
156   _ngMesh = NULL;
157 }
158
159 //================================================================================
160 /*!
161  * Set pointer to NETGENPlugin_Mesher* field of the holder, that will be
162  * nullified at destruction of this
163  */
164 //================================================================================
165
166 void NETGENPlugin_Mesher::SetSelfPointer( NETGENPlugin_Mesher ** ptr )
167 {
168   if ( _ptrToMe )
169     *_ptrToMe = NULL;
170
171   _ptrToMe = ptr;
172
173   if ( _ptrToMe )
174     *_ptrToMe = this;
175 }
176
177 //================================================================================
178 /*!
179  * \brief Initialize global NETGEN parameters with default values
180  */
181 //================================================================================
182
183 void NETGENPlugin_Mesher::SetDefaultParameters()
184 {
185   netgen::MeshingParameters& mparams = netgen::mparam;
186   // maximal mesh edge size
187   mparams.maxh = 0;//NETGENPlugin_Hypothesis::GetDefaultMaxSize();
188   mparams.minh = 0;
189   // minimal number of segments per edge
190   mparams.segmentsperedge = NETGENPlugin_Hypothesis::GetDefaultNbSegPerEdge();
191   // rate of growth of size between elements
192   mparams.grading = NETGENPlugin_Hypothesis::GetDefaultGrowthRate();
193   // safety factor for curvatures (elements per radius)
194   mparams.curvaturesafety = NETGENPlugin_Hypothesis::GetDefaultNbSegPerRadius();
195   // create elements of second order
196   mparams.secondorder = NETGENPlugin_Hypothesis::GetDefaultSecondOrder() ? 1 : 0;
197   // quad-dominated surface meshing
198   if (_isVolume)
199     mparams.quad = 0;
200   else
201     mparams.quad = NETGENPlugin_Hypothesis_2D::GetDefaultQuadAllowed() ? 1 : 0;
202   _fineness = NETGENPlugin_Hypothesis::GetDefaultFineness();
203 }
204
205 //=============================================================================
206 /*!
207  *
208  */
209 //=============================================================================
210 void SetLocalSize(TopoDS_Shape GeomShape, double LocalSize)
211 {
212   TopAbs_ShapeEnum GeomType = GeomShape.ShapeType();
213   if (GeomType == TopAbs_COMPOUND) {
214     for (TopoDS_Iterator it (GeomShape); it.More(); it.Next()) {
215       SetLocalSize(it.Value(), LocalSize);
216     }
217     return;
218   }
219   int key;
220   if (! ShapesWithLocalSize.Contains(GeomShape))
221     key = ShapesWithLocalSize.Add(GeomShape);
222   else
223     key = ShapesWithLocalSize.FindIndex(GeomShape);
224   if (GeomType == TopAbs_VERTEX) {
225     VertexId2LocalSize[key] = LocalSize;
226   } else if (GeomType == TopAbs_EDGE) {
227     EdgeId2LocalSize[key] = LocalSize;
228   } else if (GeomType == TopAbs_FACE) {
229     FaceId2LocalSize[key] = LocalSize;
230   }
231 }
232
233 //=============================================================================
234 /*!
235  * Pass parameters to NETGEN
236  */
237 //=============================================================================
238 void NETGENPlugin_Mesher::SetParameters(const NETGENPlugin_Hypothesis* hyp)
239 {
240   if (hyp)
241   {
242     netgen::MeshingParameters& mparams = netgen::mparam;
243     // Initialize global NETGEN parameters:
244     // maximal mesh segment size
245     mparams.maxh = hyp->GetMaxSize();
246     // maximal mesh element linear size
247     mparams.minh = hyp->GetMinSize();
248     // minimal number of segments per edge
249     mparams.segmentsperedge = hyp->GetNbSegPerEdge();
250     // rate of growth of size between elements
251     mparams.grading = hyp->GetGrowthRate();
252     // safety factor for curvatures (elements per radius)
253     mparams.curvaturesafety = hyp->GetNbSegPerRadius();
254     // create elements of second order
255     mparams.secondorder = hyp->GetSecondOrder() ? 1 : 0;
256     // quad-dominated surface meshing
257     // only triangles are allowed for volumic mesh (before realizing IMP 0021676)
258     //if (!_isVolume)
259       mparams.quad = hyp->GetQuadAllowed() ? 1 : 0;
260     _optimize = hyp->GetOptimize();
261     _fineness = hyp->GetFineness();
262     _simpleHyp = NULL;
263
264     SMESH_Gen_i* smeshGen_i = SMESH_Gen_i::GetSMESHGen();
265     CORBA::Object_var anObject = smeshGen_i->GetNS()->Resolve("/myStudyManager");
266     SALOMEDS::StudyManager_var aStudyMgr = SALOMEDS::StudyManager::_narrow(anObject);
267     SALOMEDS::Study_var myStudy = aStudyMgr->GetStudyByID(hyp->GetStudyId());
268     
269     const NETGENPlugin_Hypothesis::TLocalSize localSizes = hyp->GetLocalSizesAndEntries();
270     NETGENPlugin_Hypothesis::TLocalSize::const_iterator it = localSizes.begin();
271     for (it ; it != localSizes.end() ; it++)
272       {
273         std::string entry = (*it).first;
274         double val = (*it).second;
275         // --
276         GEOM::GEOM_Object_var aGeomObj;
277         TopoDS_Shape S = TopoDS_Shape();
278         SALOMEDS::SObject_var aSObj = myStudy->FindObjectID( entry.c_str() );
279         if (!aSObj->_is_nil()) {
280           CORBA::Object_var obj = aSObj->GetObject();
281           aGeomObj = GEOM::GEOM_Object::_narrow(obj);
282           aSObj->UnRegister();
283         }
284         if ( !aGeomObj->_is_nil() )
285           S = smeshGen_i->GeomObjectToShape( aGeomObj.in() );
286         // --
287         SetLocalSize(S, val);
288       }
289   }
290 }
291
292 //=============================================================================
293 /*!
294  * Pass simple parameters to NETGEN
295  */
296 //=============================================================================
297
298 void NETGENPlugin_Mesher::SetParameters(const NETGENPlugin_SimpleHypothesis_2D* hyp)
299 {
300   _simpleHyp = hyp;
301   if ( _simpleHyp )
302     SetDefaultParameters();
303 }
304
305 //=============================================================================
306 /*!
307  *  Link - a pair of integer numbers
308  */
309 //=============================================================================
310 struct Link
311 {
312   int n1, n2;
313   Link(int _n1, int _n2) : n1(_n1), n2(_n2) {}
314   Link() : n1(0), n2(0) {}
315 };
316
317 int HashCode(const Link& aLink, int aLimit)
318 {
319   return HashCode(aLink.n1 + aLink.n2, aLimit);
320 }
321
322 Standard_Boolean IsEqual(const Link& aLink1, const Link& aLink2)
323 {
324   return (aLink1.n1 == aLink2.n1 && aLink1.n2 == aLink2.n2 ||
325           aLink1.n1 == aLink2.n2 && aLink1.n2 == aLink2.n1);
326 }
327
328 namespace
329 {
330   //================================================================================
331   /*!
332    * \brief return id of netgen point corresponding to SMDS node
333    */
334   //================================================================================
335   typedef map< const SMDS_MeshNode*, int > TNode2IdMap;
336
337   int ngNodeId( const SMDS_MeshNode* node,
338                 netgen::Mesh&        ngMesh,
339                 TNode2IdMap&         nodeNgIdMap)
340   {
341     int newNgId = ngMesh.GetNP() + 1;
342
343     TNode2IdMap::iterator node_id = nodeNgIdMap.insert( make_pair( node, newNgId )).first;
344
345     if ( node_id->second == newNgId)
346     {
347 #if defined(DUMP_SEGMENTS) || defined(DUMP_TRIANGLES)
348       cout << "Ng " << newNgId << " - " << node;
349 #endif
350       netgen::MeshPoint p( netgen::Point<3> (node->X(), node->Y(), node->Z()) );
351       ngMesh.AddPoint( p );
352     }
353     return node_id->second;
354   }
355
356   //================================================================================
357   /*!
358    * \brief Return computed EDGEs connected to the given one
359    */
360   //================================================================================
361
362   list< TopoDS_Edge > getConnectedEdges( const TopoDS_Edge&                 edge,
363                                          const TopoDS_Face&                 face,
364                                          const set< SMESH_subMesh* > &      computedSM,
365                                          const SMESH_MesherHelper&          helper,
366                                          map< SMESH_subMesh*, set< int > >& addedEdgeSM2Faces)
367   {
368     // get ordered EDGEs
369     list< TopoDS_Edge > edges;
370     list< int > nbEdgesInWire;
371     int nbWires = SMESH_Block::GetOrderedEdges( face, edges, nbEdgesInWire);
372
373     // find <edge> within <edges>
374     list< TopoDS_Edge >::iterator eItFwd = edges.begin();
375     for ( ; eItFwd != edges.end(); ++eItFwd )
376       if ( edge.IsSame( *eItFwd ))
377         break;
378     if ( eItFwd == edges.end()) return list< TopoDS_Edge>();
379
380     if ( eItFwd->Orientation() >= TopAbs_INTERNAL )
381     {
382       // connected INTERNAL edges returned from GetOrderedEdges() are wrongly oriented
383       // so treat each INTERNAL edge separately
384       TopoDS_Edge e = *eItFwd;
385       edges.clear();
386       edges.push_back( e );
387       return edges;
388     }
389
390     // get all computed EDGEs connected to <edge>
391
392     list< TopoDS_Edge >::iterator eItBack = eItFwd, ePrev;
393     TopoDS_Vertex vCommon;
394     TopTools_MapOfShape eAdded; // map used not to add a seam edge twice to <edges>
395     eAdded.Add( edge );
396
397     // put edges before <edge> to <edges> back
398     while ( edges.begin() != eItFwd )
399       edges.splice( edges.end(), edges, edges.begin() );
400
401     // search forward
402     ePrev = eItFwd;
403     while ( ++eItFwd != edges.end() )
404     {
405       SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( *eItFwd );
406
407       bool connected = TopExp::CommonVertex( *ePrev, *eItFwd, vCommon );
408       bool computed  = sm->IsMeshComputed();
409       bool added     = addedEdgeSM2Faces[ sm ].count( helper.GetSubShapeID() );
410       bool doubled   = !eAdded.Add( *eItFwd );
411       bool orientOK  = (( ePrev ->Orientation() < TopAbs_INTERNAL ) ==
412                         ( eItFwd->Orientation() < TopAbs_INTERNAL )    );
413       if ( !connected || !computed || !orientOK || added || doubled )
414       {
415         // stop advancement; move edges from tail to head
416         while ( edges.back() != *ePrev )
417           edges.splice( edges.begin(), edges, --edges.end() );
418         break;
419       }
420       ePrev = eItFwd;
421     }
422     // search backward
423     while ( eItBack != edges.begin() )
424     {
425       ePrev = eItBack;
426       --eItBack;
427       SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( *eItBack );
428
429       bool connected = TopExp::CommonVertex( *ePrev, *eItBack, vCommon );
430       bool computed  = sm->IsMeshComputed();
431       bool added     = addedEdgeSM2Faces[ sm ].count( helper.GetSubShapeID() );
432       bool doubled   = !eAdded.Add( *eItBack );
433       bool orientOK  = (( ePrev  ->Orientation() < TopAbs_INTERNAL ) ==
434                         ( eItBack->Orientation() < TopAbs_INTERNAL )    );
435       if ( !connected || !computed || !orientOK || added || doubled)
436       {
437         // stop advancement
438         edges.erase( edges.begin(), ePrev );
439         break;
440       }
441     }
442     if ( edges.front() != edges.back() )
443     {
444       // assure that the 1st vertex is meshed
445       TopoDS_Edge eLast = edges.back();
446       while ( !SMESH_Algo::VertexNode( SMESH_MesherHelper::IthVertex( 0, edges.front()), helper.GetMeshDS())
447               &&
448               edges.front() != eLast )
449         edges.splice( edges.end(), edges, edges.begin() );
450     }
451     return edges;
452   }
453
454   //================================================================================
455   /*!
456    * \brief Make triangulation of a shape precise enough
457    */
458   //================================================================================
459
460   void updateTriangulation( const TopoDS_Shape& shape )
461   {
462     // static set< Poly_Triangulation* > updated;
463
464     // TopLoc_Location loc;
465     // TopExp_Explorer fExp( shape, TopAbs_FACE );
466     // for ( ; fExp.More(); fExp.Next() )
467     // {
468     //   Handle(Poly_Triangulation) triangulation =
469     //     BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
470     //   if ( triangulation.IsNull() ||
471     //        updated.insert( triangulation.operator->() ).second )
472     //   {
473     //     BRepTools::Clean (shape);
474         try {
475           OCC_CATCH_SIGNALS;
476           BRepMesh_IncrementalMesh e(shape, 0.01, true);
477         }
478         catch (Standard_Failure)
479         {
480         }
481   //       updated.erase( triangulation.operator->() );
482   //       triangulation = BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
483   //       updated.insert( triangulation.operator->() );
484   //     }
485   //   }
486   }
487 }
488
489 //================================================================================
490 /*!
491  * \brief Initialize netgen::OCCGeometry with OCCT shape
492  */
493 //================================================================================
494
495 void NETGENPlugin_Mesher::PrepareOCCgeometry(netgen::OCCGeometry&     occgeo,
496                                              const TopoDS_Shape&      shape,
497                                              SMESH_Mesh&              mesh,
498                                              list< SMESH_subMesh* > * meshedSM,
499                                              NETGENPlugin_Internals*  intern)
500 {
501   updateTriangulation( shape );
502
503   Bnd_Box bb;
504   BRepBndLib::Add (shape, bb);
505   double x1,y1,z1,x2,y2,z2;
506   bb.Get (x1,y1,z1,x2,y2,z2);
507   MESSAGE("shape bounding box:\n" <<
508           "(" << x1 << " " << y1 << " " << z1 << ") " <<
509           "(" << x2 << " " << y2 << " " << z2 << ")");
510   netgen::Point<3> p1 = netgen::Point<3> (x1,y1,z1);
511   netgen::Point<3> p2 = netgen::Point<3> (x2,y2,z2);
512   occgeo.boundingbox = netgen::Box<3> (p1,p2);
513
514   occgeo.shape = shape;
515   occgeo.changed = 1;
516
517   // fill maps of shapes of occgeo with not yet meshed subshapes
518
519   // get root submeshes
520   list< SMESH_subMesh* > rootSM;
521   const int shapeID = mesh.GetMeshDS()->ShapeToIndex( shape );
522   if ( shapeID > 0 ) { // SMESH_subMesh with ID 0 may exist, don't use it!
523     rootSM.push_back( mesh.GetSubMesh( shape ));
524   }
525   else {
526     for ( TopoDS_Iterator it( shape ); it.More(); it.Next() )
527       rootSM.push_back( mesh.GetSubMesh( it.Value() ));
528   }
529
530   // add subshapes of empty submeshes
531   list< SMESH_subMesh* >::iterator rootIt = rootSM.begin(), rootEnd = rootSM.end();
532   for ( ; rootIt != rootEnd; ++rootIt ) {
533     SMESH_subMesh * root = *rootIt;
534     SMESH_subMeshIteratorPtr smIt = root->getDependsOnIterator(/*includeSelf=*/true,
535                                                                /*complexShapeFirst=*/true);
536     // to find a right orientation of subshapes (PAL20462)
537     TopTools_IndexedMapOfShape subShapes;
538     TopExp::MapShapes(root->GetSubShape(), subShapes);
539     while ( smIt->more() )
540     {
541       SMESH_subMesh* sm = smIt->next();
542       TopoDS_Shape shape = sm->GetSubShape();
543       if ( intern && intern->isShapeToPrecompute( shape ))
544         continue;
545       if ( !meshedSM || sm->IsEmpty() )
546       {
547         if ( shape.ShapeType() != TopAbs_VERTEX )
548           shape = subShapes( subShapes.FindIndex( shape ));// shape -> index -> oriented shape
549         if ( shape.Orientation() >= TopAbs_INTERNAL )
550           shape.Orientation( TopAbs_FORWARD ); // isuue 0020676
551         switch ( shape.ShapeType() ) {
552         case TopAbs_FACE  : occgeo.fmap.Add( shape ); break;
553         case TopAbs_EDGE  : occgeo.emap.Add( shape ); break;
554         case TopAbs_VERTEX: occgeo.vmap.Add( shape ); break;
555         case TopAbs_SOLID :occgeo.somap.Add( shape ); break;
556         default:;
557         }
558       }
559       // collect submeshes of meshed shapes
560       else if (meshedSM)
561       {
562         const int dim = SMESH_Gen::GetShapeDim( shape );
563         meshedSM[ dim ].push_back( sm );
564       }
565     }
566   }
567   occgeo.facemeshstatus.SetSize (occgeo.fmap.Extent());
568   occgeo.facemeshstatus = 0;
569   occgeo.face_maxh_modified.SetSize(occgeo.fmap.Extent());
570   occgeo.face_maxh_modified = 0;
571   occgeo.face_maxh.SetSize(occgeo.fmap.Extent());
572   occgeo.face_maxh = netgen::mparam.maxh;
573 }
574
575 //================================================================================
576 /*!
577  * \brief Return a default min size value suitable for the given geometry.
578  */
579 //================================================================================
580
581 double NETGENPlugin_Mesher::GetDefaultMinSize(const TopoDS_Shape& geom,
582                                               const double        maxSize)
583 {
584   updateTriangulation( geom );
585
586   TopLoc_Location loc;
587   int i1, i2, i3;
588   const int* pi[4] = { &i1, &i2, &i3, &i1 };
589   double minh = 1e100;
590   Bnd_B3d bb;
591   TopExp_Explorer fExp( geom, TopAbs_FACE );
592   for ( ; fExp.More(); fExp.Next() )
593   {
594     Handle(Poly_Triangulation) triangulation =
595       BRep_Tool::Triangulation ( TopoDS::Face( fExp.Current() ), loc);
596     if ( triangulation.IsNull() ) continue;
597     const double fTol = BRep_Tool::Tolerance( TopoDS::Face( fExp.Current() ));
598     const TColgp_Array1OfPnt&   points = triangulation->Nodes();
599     const Poly_Array1OfTriangle& trias = triangulation->Triangles();
600     for ( int iT = trias.Lower(); iT <= trias.Upper(); ++iT )
601     {
602       trias(iT).Get( i1, i2, i3 );
603       for ( int j = 0; j < 3; ++j )
604       {
605         double dist2 = points(*pi[j]).SquareDistance( points( *pi[j+1] ));
606         if ( dist2 < minh && fTol*fTol < dist2 )
607           minh = dist2;
608         bb.Add( points(*pi[j]));
609       }
610     }
611   }
612   if ( minh > 0.25 * bb.SquareExtent() ) // simple geometry, rough triangulation
613   {
614     minh = 1e-3 * sqrt( bb.SquareExtent());
615     //cout << "BND BOX minh = " <<minh << endl;
616   }
617   else
618   {
619     minh = 3 * sqrt( minh ); // triangulation for visualization is rather fine
620     //cout << "TRIANGULATION minh = " <<minh << endl;
621   }
622   if ( minh > 0.5 * maxSize )
623     minh = maxSize / 3.;
624
625   return minh;
626 }
627
628 //================================================================================
629 /*!
630  * \brief Restrict size of elements at a given point
631  */
632 //================================================================================
633
634 void NETGENPlugin_Mesher::RestrictLocalSize(netgen::Mesh& ngMesh, const gp_XYZ& p, const double size)
635 {
636   if ( netgen::mparam.minh > size )
637   {
638     ngMesh.SetMinimalH( size );
639     netgen::mparam.minh = size;
640   }
641   netgen::Point3d pi(p.X(), p.Y(), p.Z());
642   ngMesh.RestrictLocalH( pi, size );
643 }
644
645 //================================================================================
646 /*!
647  * \brief fill ngMesh with nodes and elements of computed submeshes
648  */
649 //================================================================================
650
651 bool NETGENPlugin_Mesher::FillNgMesh(netgen::OCCGeometry&           occgeom,
652                                      netgen::Mesh&                  ngMesh,
653                                      vector<const SMDS_MeshNode*>&  nodeVec,
654                                      const list< SMESH_subMesh* > & meshedSM,
655                                      SMESH_ProxyMesh::Ptr           proxyMesh)
656 {
657   TNode2IdMap nodeNgIdMap;
658   for ( int i = 1; i < nodeVec.size(); ++i )
659     nodeNgIdMap.insert( make_pair( nodeVec[i], i ));
660
661   TopTools_MapOfShape visitedShapes;
662   map< SMESH_subMesh*, set< int > > visitedEdgeSM2Faces;
663   set< SMESH_subMesh* > computedSM( meshedSM.begin(), meshedSM.end() );
664
665   SMESH_MesherHelper helper (*_mesh);
666
667   int faceNgID = ngMesh.GetNFD();
668
669   list< SMESH_subMesh* >::const_iterator smIt, smEnd = meshedSM.end();
670   for ( smIt = meshedSM.begin(); smIt != smEnd; ++smIt )
671   {
672     SMESH_subMesh* sm = *smIt;
673     if ( !visitedShapes.Add( sm->GetSubShape() ))
674       continue;
675
676     const SMESHDS_SubMesh * smDS = sm->GetSubMeshDS();
677     if ( !smDS ) continue;
678
679     switch ( sm->GetSubShape().ShapeType() )
680     {
681     case TopAbs_EDGE: { // EDGE
682       // ----------------------
683       TopoDS_Edge geomEdge  = TopoDS::Edge( sm->GetSubShape() );
684       if ( geomEdge.Orientation() >= TopAbs_INTERNAL )
685         geomEdge.Orientation( TopAbs_FORWARD ); // issue 0020676
686
687       // Add ng segments for each not meshed FACE the EDGE bounds
688       PShapeIteratorPtr fIt = helper.GetAncestors( geomEdge, *sm->GetFather(), TopAbs_FACE );
689       while ( const TopoDS_Shape * anc = fIt->next() )
690       {
691         faceNgID = occgeom.fmap.FindIndex( *anc );
692         if ( faceNgID < 1 )
693           continue; // meshed face
694
695         int faceSMDSId = helper.GetMeshDS()->ShapeToIndex( *anc );
696         if ( visitedEdgeSM2Faces[ sm ].count( faceSMDSId ))
697           continue; // already treated EDGE
698
699         TopoDS_Face face = TopoDS::Face( occgeom.fmap( faceNgID ));
700         if ( face.Orientation() >= TopAbs_INTERNAL )
701           face.Orientation( TopAbs_FORWARD ); // issue 0020676
702
703         // get all meshed EDGEs of the FACE connected to geomEdge (issue 0021140)
704         helper.SetSubShape( face );
705         list< TopoDS_Edge > edges = getConnectedEdges( geomEdge, face, computedSM, helper,
706                                                        visitedEdgeSM2Faces );
707         if ( edges.empty() )
708           continue; // wrong ancestor?
709
710         // find out orientation of <edges> within <face>
711         TopoDS_Edge eNotSeam = edges.front();
712         if ( helper.HasSeam() )
713         {
714           list< TopoDS_Edge >::iterator eIt = edges.begin();
715           while ( helper.IsRealSeam( *eIt )) ++eIt;
716           if ( eIt != edges.end() )
717             eNotSeam = *eIt;
718         }
719         TopAbs_Orientation fOri = helper.GetSubShapeOri( face, eNotSeam );
720         bool isForwad = ( fOri == eNotSeam.Orientation() || fOri >= TopAbs_INTERNAL );
721
722         // get all nodes from connected <edges>
723         bool isQuad   = smDS->NbElements() ? smDS->GetElements()->next()->IsQuadratic() : false;
724         StdMeshers_FaceSide fSide( face, edges, _mesh, isForwad, isQuad );
725         const vector<UVPtStruct>& points = fSide.GetUVPtStruct();
726         int i, nbSeg = fSide.NbSegments();
727
728         // remember EDGEs of fSide to treat only once
729         for ( int iE = 0; iE < fSide.NbEdges(); ++iE )
730           visitedEdgeSM2Faces[ helper.GetMesh()->GetSubMesh( fSide.Edge(iE )) ].insert(faceSMDSId);
731
732         double otherSeamParam = 0;
733         bool isSeam = false;
734
735         // add segments
736
737         int prevNgId = ngNodeId( points[0].node, ngMesh, nodeNgIdMap );
738
739         for ( i = 0; i < nbSeg; ++i )
740         {
741           const UVPtStruct& p1 = points[ i ];
742           const UVPtStruct& p2 = points[ i+1 ];
743
744           if ( p1.node->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX ) //an EDGE begins
745           {
746             isSeam = false;
747             if ( helper.IsRealSeam( p1.node->getshapeId() ))
748             {
749               TopoDS_Edge e = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
750               isSeam = helper.IsRealSeam( e );
751               if ( isSeam )
752               {
753                 otherSeamParam = helper.GetOtherParam( helper.GetPeriodicIndex() & 1 ? p2.u : p2.v );
754               }
755             }
756           }
757           netgen::Segment seg;
758           // ng node ids
759           seg[0] = prevNgId;
760           seg[1] = prevNgId = ngNodeId( p2.node, ngMesh, nodeNgIdMap );
761           // node param on curve
762           seg.epgeominfo[ 0 ].dist = p1.param;
763           seg.epgeominfo[ 1 ].dist = p2.param;
764           // uv on face
765           seg.epgeominfo[ 0 ].u = p1.u;
766           seg.epgeominfo[ 0 ].v = p1.v;
767           seg.epgeominfo[ 1 ].u = p2.u;
768           seg.epgeominfo[ 1 ].v = p2.v;
769
770           //geomEdge = fSide.Edge( fSide.EdgeIndex( 0.5 * ( p1.normParam + p2.normParam )));
771           //seg.epgeominfo[ 0 ].edgenr = seg.epgeominfo[ 1 ].edgenr = occgeom.emap.FindIndex( geomEdge );
772
773           //seg.epgeominfo[ iEnd ].edgenr = edgeID; //  = geom.emap.FindIndex(edge);
774           seg.si = faceNgID;                   // = geom.fmap.FindIndex (face);
775           seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
776           ngMesh.AddSegment (seg);
777
778           SMESH_TNodeXYZ np1( p1.node ), np2( p2.node );
779           RestrictLocalSize( ngMesh, 0.5*(np1+np2), (np1-np2).Modulus() );
780
781 #ifdef DUMP_SEGMENTS
782           cout << "Segment: " << seg.edgenr << " on SMESH face " << helper.GetMeshDS()->ShapeToIndex( face ) << endl
783                << "\tface index: " << seg.si << endl
784                << "\tp1: " << seg[0] << endl
785                << "\tp2: " << seg[1] << endl
786                << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
787                << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
788             //<< "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
789                << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
790                << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl;
791             //<< "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
792 #endif
793           if ( isSeam )
794           {
795             if ( helper.GetPeriodicIndex() && 1 ) {
796               seg.epgeominfo[ 0 ].u = otherSeamParam;
797               seg.epgeominfo[ 1 ].u = otherSeamParam;
798               swap (seg.epgeominfo[0].v, seg.epgeominfo[1].v);
799             } else {
800               seg.epgeominfo[ 0 ].v = otherSeamParam;
801               seg.epgeominfo[ 1 ].v = otherSeamParam;
802               swap (seg.epgeominfo[0].u, seg.epgeominfo[1].u);
803             }
804             swap (seg[0], seg[1]);
805             swap (seg.epgeominfo[0].dist, seg.epgeominfo[1].dist);
806             seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
807             ngMesh.AddSegment (seg);
808 #ifdef DUMP_SEGMENTS
809             cout << "Segment: " << seg.edgenr << endl
810                  << "\t is SEAM (reverse) of the previous. "
811                  << " Other " << (helper.GetPeriodicIndex() && 1 ? "U" : "V")
812                  << " = " << otherSeamParam << endl;
813 #endif
814           }
815           else if ( fOri == TopAbs_INTERNAL )
816           {
817             swap (seg[0], seg[1]);
818             swap( seg.epgeominfo[0], seg.epgeominfo[1] );
819             seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
820             ngMesh.AddSegment (seg);
821 #ifdef DUMP_SEGMENTS
822             cout << "Segment: " << seg.edgenr << endl << "\t is REVERSE of the previous" << endl;
823 #endif
824           }
825         }
826       } // loop on geomEdge ancestors
827
828       break;
829     } // case TopAbs_EDGE
830
831     case TopAbs_FACE: { // FACE
832       // ----------------------
833       const TopoDS_Face& geomFace  = TopoDS::Face( sm->GetSubShape() );
834       helper.SetSubShape( geomFace );
835       bool isInternalFace = ( geomFace.Orientation() == TopAbs_INTERNAL );
836
837       // Find solids the geomFace bounds
838       int solidID1 = 0, solidID2 = 0;
839       StdMeshers_QuadToTriaAdaptor* quadAdaptor =
840         dynamic_cast<StdMeshers_QuadToTriaAdaptor*>( proxyMesh.get() );
841       if ( quadAdaptor )
842       {
843         solidID1 = occgeom.somap.FindIndex( quadAdaptor->GetShape() );
844       }
845       else
846       {  
847         PShapeIteratorPtr solidIt = helper.GetAncestors( geomFace, *sm->GetFather(), TopAbs_SOLID);
848         while ( const TopoDS_Shape * solid = solidIt->next() )
849         {
850           int id = occgeom.somap.FindIndex ( *solid );
851           if ( solidID1 && id != solidID1 ) solidID2 = id;
852           else                              solidID1 = id;
853         }
854       }
855       // Add ng face descriptors of meshed faces
856       faceNgID++;
857       ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(faceNgID, solidID1, solidID2, 0));
858
859       // if second oreder is required, even already meshed faces must be passed to NETGEN
860       int fID = occgeom.fmap.Add( geomFace );
861       while ( fID < faceNgID ) // geomFace is already in occgeom.fmap, add a copy
862         fID = occgeom.fmap.Add( BRepBuilderAPI_Copy( geomFace, /*copyGeom=*/false ));
863       // Problem with the second order in a quadrangular mesh remains.
864       // 1) All quadrangles geberated by NETGEN are moved to an inexistent face
865       //    by FillSMesh() (find AddFaceDescriptor)
866       // 2) Temporary triangles generated by StdMeshers_QuadToTriaAdaptor
867       //    are on faces where quadrangles were.
868       // Due to these 2 points, wrong geom faces are used while conversion to qudratic
869       // of the mentioned above quadrangles and triangles
870
871       // Orient the face correctly in solidID1 (issue 0020206)
872       bool reverse = false;
873       if ( solidID1 ) {
874         TopoDS_Shape solid = occgeom.somap( solidID1 );
875         TopAbs_Orientation faceOriInSolid = helper.GetSubShapeOri( solid, geomFace );
876         if ( faceOriInSolid >= 0 )
877           reverse =
878             helper.IsReversedSubMesh( TopoDS::Face( geomFace.Oriented( faceOriInSolid )));
879       }
880
881       // Add surface elements
882
883       netgen::Element2d tri(3);
884       tri.SetIndex ( faceNgID );
885
886
887 #ifdef DUMP_TRIANGLES
888       cout << "SMESH face " << helper.GetMeshDS()->ShapeToIndex( geomFace )
889            << " internal="<<isInternalFace << endl;
890 #endif
891       if ( proxyMesh )
892         smDS = proxyMesh->GetSubMesh( geomFace );
893
894       SMDS_ElemIteratorPtr faces = smDS->GetElements();
895       while ( faces->more() )
896       {
897         const SMDS_MeshElement* f = faces->next();
898         if ( f->NbNodes() % 3 != 0 ) // not triangle
899         {
900           PShapeIteratorPtr solidIt=helper.GetAncestors(geomFace,*sm->GetFather(),TopAbs_SOLID);
901           if ( const TopoDS_Shape * solid = solidIt->next() )
902             sm = _mesh->GetSubMesh( *solid );
903           SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
904           smError.reset( new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH,"Not triangle submesh"));
905           smError->myBadElements.push_back( f );
906           return false;
907         }
908
909         for ( int i = 0; i < 3; ++i )
910         {
911           const SMDS_MeshNode* node = f->GetNode( i ), * inFaceNode=0;
912
913           // get node UV on face
914           int shapeID = node->getshapeId();
915           if ( helper.IsSeamShape( shapeID ))
916             if ( helper.IsSeamShape( f->GetNodeWrap( i+1 )->getshapeId() ))
917               inFaceNode = f->GetNodeWrap( i-1 );
918             else 
919               inFaceNode = f->GetNodeWrap( i+1 );
920           gp_XY uv = helper.GetNodeUV( geomFace, node, inFaceNode );
921
922           int ind = reverse ? 3-i : i+1;
923           tri.GeomInfoPi(ind).u = uv.X();
924           tri.GeomInfoPi(ind).v = uv.Y();
925           tri.PNum      (ind) = ngNodeId( node, ngMesh, nodeNgIdMap );
926         }
927
928         ngMesh.AddSurfaceElement (tri);
929 #ifdef DUMP_TRIANGLES
930         cout << tri << endl;
931 #endif
932
933         if ( isInternalFace )
934         {
935           swap( tri[1], tri[2] );
936           ngMesh.AddSurfaceElement (tri);
937 #ifdef DUMP_TRIANGLES
938         cout << tri << endl;
939 #endif
940         }
941       }
942       break;
943     } // case TopAbs_FACE
944
945     case TopAbs_VERTEX: { // VERTEX
946       // --------------------------
947       // issue 0021405. Add node only if a VERTEX is shared by a not meshed EDGE,
948       // else netgen removes a free node and nodeVector becomes invalid
949       PShapeIteratorPtr ansIt = helper.GetAncestors( sm->GetSubShape(),
950                                                      *sm->GetFather(),
951                                                      TopAbs_EDGE );
952       bool toAdd = false;
953       while ( const TopoDS_Shape* e = ansIt->next() )
954       {
955         SMESH_subMesh* eSub = helper.GetMesh()->GetSubMesh( *e );
956         if (( toAdd = eSub->IsEmpty() )) break;
957       }
958       if ( toAdd )
959       {
960         SMDS_NodeIteratorPtr nodeIt = smDS->GetNodes();
961         if ( nodeIt->more() )
962           ngNodeId( nodeIt->next(), ngMesh, nodeNgIdMap );
963       }
964       break;
965     }
966     default:;
967     } // switch
968   } // loop on submeshes
969
970   // fill nodeVec
971   nodeVec.resize( ngMesh.GetNP() + 1 );
972   TNode2IdMap::iterator node_NgId, nodeNgIdEnd = nodeNgIdMap.end();
973   for ( node_NgId = nodeNgIdMap.begin(); node_NgId != nodeNgIdEnd; ++node_NgId)
974     nodeVec[ node_NgId->second ] = node_NgId->first;
975
976   return true;
977 }
978
979 //================================================================================
980 /*!
981  * \brief Duplicate mesh faces on internal geom faces
982  */
983 //================================================================================
984
985 void NETGENPlugin_Mesher::FixIntFaces(const netgen::OCCGeometry& occgeom,
986                                       netgen::Mesh&              ngMesh,
987                                       NETGENPlugin_Internals&    internalShapes)
988 {
989   SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
990   
991   // find ng indices of internal faces
992   set<int> ngFaceIds;
993   for ( int ngFaceID = 1; ngFaceID <= occgeom.fmap.Extent(); ++ngFaceID )
994   {
995     int smeshID = meshDS->ShapeToIndex( occgeom.fmap( ngFaceID ));
996     if ( internalShapes.isInternalShape( smeshID ))
997       ngFaceIds.insert( ngFaceID );
998   }
999   if ( !ngFaceIds.empty() )
1000   {
1001     // duplicate faces
1002     int i, nbFaces = ngMesh.GetNSE();
1003     for (int i = 1; i <= nbFaces; ++i)
1004     {
1005       netgen::Element2d elem = ngMesh.SurfaceElement(i);
1006       if ( ngFaceIds.count( elem.GetIndex() ))
1007       {
1008         swap( elem[1], elem[2] );
1009         ngMesh.AddSurfaceElement (elem);
1010       }
1011     }
1012   }
1013 }
1014
1015 namespace
1016 {
1017   //================================================================================
1018   // define gp_XY_Subtracted pointer to function calling gp_XY::Subtracted(gp_XY)
1019   gp_XY_FunPtr(Subtracted);
1020   //gp_XY_FunPtr(Added);
1021
1022   //================================================================================
1023   /*!
1024    * \brief Evaluate distance between two 2d points along the surface
1025    */
1026   //================================================================================
1027
1028   double evalDist( const gp_XY&                uv1,
1029                    const gp_XY&                uv2,
1030                    const Handle(Geom_Surface)& surf,
1031                    const int                   stopHandler=-1)
1032   {
1033     if ( stopHandler > 0 ) // continue recursion
1034     {
1035       gp_XY mid = SMESH_MesherHelper::GetMiddleUV( surf, uv1, uv2 );
1036       return evalDist( uv1,mid, surf, stopHandler-1 ) + evalDist( mid,uv2, surf, stopHandler-1 );
1037     }
1038     double dist3D = surf->Value( uv1.X(), uv1.Y() ).Distance( surf->Value( uv2.X(), uv2.Y() ));
1039     if ( stopHandler == 0 ) // stop recursion
1040       return dist3D;
1041     
1042     // start recursion if necessary
1043     double dist2D = SMESH_MesherHelper::applyIn2D(surf, uv1, uv2, gp_XY_Subtracted, 0).Modulus();
1044     if ( fabs( dist3D - dist2D ) < dist2D * 1e-10 )
1045       return dist3D; // equal parametrization of a planar surface
1046
1047     return evalDist( uv1, uv2, surf, 3 ); // start recursion
1048   }
1049
1050   //================================================================================
1051   /*!
1052    * \brief Data of vertex internal in geom face
1053    */
1054   //================================================================================
1055
1056   struct TIntVData
1057   {
1058     gp_XY uv;        //!< UV in face parametric space
1059     int   ngId;      //!< ng id of corrsponding node
1060     gp_XY uvClose;   //!< UV of closest boundary node
1061     int   ngIdClose; //!< ng id of closest boundary node
1062   };
1063
1064   //================================================================================
1065   /*!
1066    * \brief Data of vertex internal in solid
1067    */
1068   //================================================================================
1069
1070   struct TIntVSoData
1071   {
1072     int   ngId;      //!< ng id of corresponding node
1073     int   ngIdClose; //!< ng id of closest 2d mesh element
1074     int   ngIdCloseN; //!< ng id of closest node of the closest 2d mesh element
1075   };
1076
1077   inline double dist2(const netgen::MeshPoint& p1, const netgen::MeshPoint& p2)
1078   {
1079     return gp_Pnt( NGPOINT_COORDS(p1)).SquareDistance( gp_Pnt( NGPOINT_COORDS(p2)));
1080   }
1081 }
1082
1083 //================================================================================
1084 /*!
1085  * \brief Make netgen take internal vertices in faces into account by adding
1086  *        segments including internal vertices
1087  *
1088  * This function works in supposition that 1D mesh is already computed in ngMesh
1089  */
1090 //================================================================================
1091
1092 void NETGENPlugin_Mesher::AddIntVerticesInFaces(const netgen::OCCGeometry&     occgeom,
1093                                                 netgen::Mesh&                  ngMesh,
1094                                                 vector<const SMDS_MeshNode*>&  nodeVec,
1095                                                 NETGENPlugin_Internals&        internalShapes)
1096 {
1097   if ( nodeVec.size() < ngMesh.GetNP() )
1098     nodeVec.resize( ngMesh.GetNP(), 0 );
1099
1100   SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1101   SMESH_MesherHelper helper( internalShapes.getMesh() );
1102
1103   const map<int,list<int> >& face2Vert = internalShapes.getFacesWithVertices();
1104   map<int,list<int> >::const_iterator f2v = face2Vert.begin();
1105   for ( ; f2v != face2Vert.end(); ++f2v )
1106   {
1107     const TopoDS_Face& face = TopoDS::Face( meshDS->IndexToShape( f2v->first ));
1108     if ( face.IsNull() ) continue;
1109     int faceNgID = occgeom.fmap.FindIndex (face);
1110     if ( faceNgID < 0 ) continue;
1111
1112     TopLoc_Location loc;
1113     Handle(Geom_Surface) surf = BRep_Tool::Surface(face,loc);
1114
1115     helper.SetSubShape( face );
1116     helper.SetElementsOnShape( true );
1117
1118     // Get data of internal vertices and add them to ngMesh
1119
1120     multimap< double, TIntVData > dist2VData; // sort vertices by distance from boundary nodes
1121
1122     int i, nbSegInit = ngMesh.GetNSeg();
1123
1124     // boundary characteristics
1125     double totSegLen2D = 0;
1126     int totNbSeg = 0;
1127
1128     const list<int>& iVertices = f2v->second;
1129     list<int>::const_iterator iv = iVertices.begin();
1130     for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1131     {
1132       TIntVData vData;
1133       // get node on vertex
1134       const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1135       const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1136       if ( !nV )
1137       {
1138         SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1139         sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1140         nV = SMESH_Algo::VertexNode( V, meshDS );
1141         if ( !nV ) continue;
1142       }
1143       // add ng node
1144       netgen::MeshPoint mp( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1145       ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1146       vData.ngId = ngMesh.GetNP();
1147       nodeVec.push_back( nV );
1148
1149       // get node UV
1150       bool uvOK = false;
1151       vData.uv = helper.GetNodeUV( face, nV, 0, &uvOK );
1152       if ( !uvOK ) helper.CheckNodeUV( face, nV, vData.uv, BRep_Tool::Tolerance(V),/*force=*/1);
1153
1154       // loop on all segments of the face to find the node closest to vertex and to count
1155       // average segment 2d length
1156       double closeDist2 = numeric_limits<double>::max(), dist2;
1157       int ngIdLast = 0;
1158       for (i = 1; i <= ngMesh.GetNSeg(); ++i)
1159       {
1160         netgen::Segment & seg = ngMesh.LineSegment(i);
1161         if ( seg.si != faceNgID ) continue;
1162         gp_XY uv[2];
1163         for ( int iEnd = 0; iEnd < 2; ++iEnd)
1164         {
1165           uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1166           if ( ngIdLast == seg[ iEnd ] ) continue;
1167           dist2 = helper.applyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1168           if ( dist2 < closeDist2 )
1169             vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1170           ngIdLast = seg[ iEnd ];
1171         }
1172         if ( !nbV )
1173         {
1174           totSegLen2D += helper.applyIn2D(surf, uv[0], uv[1], gp_XY_Subtracted, false).Modulus();
1175           totNbSeg++;
1176         }
1177       }
1178       dist2VData.insert( make_pair( closeDist2, vData ));
1179     }
1180
1181     if ( totNbSeg == 0 ) break;
1182     double avgSegLen2d = totSegLen2D / totNbSeg;
1183
1184     // Loop on vertices to add segments
1185
1186     multimap< double, TIntVData >::iterator dist_vData = dist2VData.begin();
1187     for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1188     {
1189       double closeDist2 = dist_vData->first, dist2;
1190       TIntVData & vData = dist_vData->second;
1191
1192       // try to find more close node among segments added for internal vertices
1193       for (i = nbSegInit+1; i <= ngMesh.GetNSeg(); ++i)
1194       {
1195         netgen::Segment & seg = ngMesh.LineSegment(i);
1196         if ( seg.si != faceNgID ) continue;
1197         gp_XY uv[2];
1198         for ( int iEnd = 0; iEnd < 2; ++iEnd)
1199         {
1200           uv[iEnd].SetCoord( seg.epgeominfo[iEnd].u, seg.epgeominfo[iEnd].v );
1201           dist2 = helper.applyIn2D(surf, uv[iEnd], vData.uv, gp_XY_Subtracted,0).SquareModulus();
1202           if ( dist2 < closeDist2 )
1203             vData.ngIdClose = seg[ iEnd ], vData.uvClose = uv[iEnd], closeDist2 = dist2;
1204         }
1205       }
1206       // decide whether to use the closest node as the second end of segment or to
1207       // create a new point
1208       int segEnd1 = vData.ngId;
1209       int segEnd2 = vData.ngIdClose; // to use closest node
1210       gp_XY uvV = vData.uv, uvP = vData.uvClose;
1211       double segLenHint  = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1212       double nodeDist2D  = sqrt( closeDist2 );
1213       double nodeDist3D  = evalDist( vData.uv, vData.uvClose, surf );
1214       bool avgLenOK  = ( avgSegLen2d < 0.75 * nodeDist2D );
1215       bool hintLenOK = ( segLenHint  < 0.75 * nodeDist3D );
1216       //cout << "uvV " << uvV.X() <<","<<uvV.Y() << " ";
1217       if ( hintLenOK || avgLenOK )
1218       {
1219         // create a point between the closest node and V
1220
1221         // how far from V
1222         double r = min( 0.5, ( hintLenOK ? segLenHint/nodeDist3D : avgSegLen2d/nodeDist2D ));
1223         // direction from V to closet node in 2D
1224         gp_Dir2d v2n( helper.applyIn2D(surf, uvP, uvV, gp_XY_Subtracted, false ));
1225         // new point
1226         uvP = vData.uv + r * nodeDist2D * v2n.XY();
1227         gp_Pnt P = surf->Value( uvP.X(), uvP.Y() ).Transformed( loc );
1228
1229         netgen::MeshPoint mp( netgen::Point<3> (P.X(), P.Y(), P.Z()));
1230         ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1231         segEnd2 = ngMesh.GetNP();
1232         //cout << "Middle " << r << " uv " << uvP.X() << "," << uvP.Y() << "( " << ngMesh.Point(segEnd2).X()<<","<<ngMesh.Point(segEnd2).Y()<<","<<ngMesh.Point(segEnd2).Z()<<" )"<< endl;
1233         SMDS_MeshNode * nP = helper.AddNode(P.X(), P.Y(), P.Z());
1234         nodeVec.push_back( nP );
1235       }
1236       //else cout << "at Node " << " uv " << uvP.X() << "," << uvP.Y() << endl;
1237
1238       // Add the segment
1239       netgen::Segment seg;
1240
1241       if ( segEnd1 > segEnd2 ) swap( segEnd1, segEnd2 ), swap( uvV, uvP );
1242       seg[0] = segEnd1;  // ng node id
1243       seg[1] = segEnd2;  // ng node id
1244       seg.edgenr = ngMesh.GetNSeg() + 1;// segment id
1245       seg.si = faceNgID;
1246
1247       seg.epgeominfo[ 0 ].dist = 0; // param on curve
1248       seg.epgeominfo[ 0 ].u    = uvV.X();
1249       seg.epgeominfo[ 0 ].v    = uvV.Y();
1250       seg.epgeominfo[ 1 ].dist = 1; // param on curve
1251       seg.epgeominfo[ 1 ].u    = uvP.X();
1252       seg.epgeominfo[ 1 ].v    = uvP.Y();
1253
1254 //       seg.epgeominfo[ 0 ].edgenr = 10; //  = geom.emap.FindIndex(edge);
1255 //       seg.epgeominfo[ 1 ].edgenr = 10; //  = geom.emap.FindIndex(edge);
1256
1257       ngMesh.AddSegment (seg);
1258
1259       // add reverse segment
1260       swap (seg[0], seg[1]);
1261       swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1262       seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1263       ngMesh.AddSegment (seg);
1264     }
1265
1266   }
1267 }
1268
1269 //================================================================================
1270 /*!
1271  * \brief Make netgen take internal vertices in solids into account by adding
1272  *        faces including internal vertices
1273  *
1274  * This function works in supposition that 2D mesh is already computed in ngMesh
1275  */
1276 //================================================================================
1277
1278 void NETGENPlugin_Mesher::AddIntVerticesInSolids(const netgen::OCCGeometry&     occgeom,
1279                                                  netgen::Mesh&                  ngMesh,
1280                                                  vector<const SMDS_MeshNode*>&  nodeVec,
1281                                                  NETGENPlugin_Internals&        internalShapes)
1282 {
1283 #ifdef DUMP_TRIANGLES_SCRIPT
1284   // create a python script making a mesh containing triangles added for internal vertices
1285   ofstream py(DUMP_TRIANGLES_SCRIPT);
1286   py << "from smesh import * "<< endl
1287      << "m = Mesh(name='triangles')" << endl;
1288 #endif
1289   if ( nodeVec.size() < ngMesh.GetNP() )
1290     nodeVec.resize( ngMesh.GetNP(), 0 );
1291
1292   SMESHDS_Mesh* meshDS = internalShapes.getMesh().GetMeshDS();
1293   SMESH_MesherHelper helper( internalShapes.getMesh() );
1294
1295   const map<int,list<int> >& so2Vert = internalShapes.getSolidsWithVertices();
1296   map<int,list<int> >::const_iterator s2v = so2Vert.begin();
1297   for ( ; s2v != so2Vert.end(); ++s2v )
1298   {
1299     const TopoDS_Shape& solid = meshDS->IndexToShape( s2v->first );
1300     if ( solid.IsNull() ) continue;
1301     int solidNgID = occgeom.somap.FindIndex (solid);
1302     if ( solidNgID < 0 && !occgeom.somap.IsEmpty() ) continue;
1303
1304     helper.SetSubShape( solid );
1305     helper.SetElementsOnShape( true );
1306
1307     // find ng indices of faces within the solid
1308     set<int> ngFaceIds;
1309     for (TopExp_Explorer fExp(solid, TopAbs_FACE); fExp.More(); fExp.Next() )
1310       ngFaceIds.insert( occgeom.fmap.FindIndex( fExp.Current() ));
1311     if ( ngFaceIds.size() == 1 && *ngFaceIds.begin() == 0 )
1312       ngFaceIds.insert( 1 );
1313
1314     // Get data of internal vertices and add them to ngMesh
1315
1316     multimap< double, TIntVSoData > dist2VData; // sort vertices by distance from ng faces
1317
1318     int i, nbFaceInit = ngMesh.GetNSE();
1319
1320     // boundary characteristics
1321     double totSegLen = 0;
1322     int totNbSeg = 0;
1323
1324     const list<int>& iVertices = s2v->second;
1325     list<int>::const_iterator iv = iVertices.begin();
1326     for ( int nbV = 0; iv != iVertices.end(); ++iv, nbV++ )
1327     {
1328       TIntVSoData vData;
1329       const TopoDS_Vertex V = TopoDS::Vertex( meshDS->IndexToShape( *iv ));
1330
1331       // get node on vertex
1332       const SMDS_MeshNode * nV = SMESH_Algo::VertexNode( V, meshDS );
1333       if ( !nV )
1334       {
1335         SMESH_subMesh* sm = helper.GetMesh()->GetSubMesh( V );
1336         sm->ComputeStateEngine( SMESH_subMesh::COMPUTE );
1337         nV = SMESH_Algo::VertexNode( V, meshDS );
1338         if ( !nV ) continue;
1339       }
1340       // add ng node
1341       netgen::MeshPoint mpV( netgen::Point<3> (nV->X(), nV->Y(), nV->Z()) );
1342       ngMesh.AddPoint ( mpV, 1, netgen::FIXEDPOINT );
1343       vData.ngId = ngMesh.GetNP();
1344       nodeVec.push_back( nV );
1345
1346       // loop on all 2d elements to find the one closest to vertex and to count
1347       // average segment length
1348       double closeDist2 = numeric_limits<double>::max(), avgDist2;
1349       for (i = 1; i <= ngMesh.GetNSE(); ++i)
1350       {
1351         const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1352         if ( !ngFaceIds.count( elem.GetIndex() )) continue;
1353         avgDist2 = 0;
1354         multimap< double, int> dist2nID; // sort nodes of element by distance from V
1355         for ( int j = 0; j < elem.GetNP(); ++j)
1356         {
1357           netgen::MeshPoint mp = ngMesh.Point( elem[j] );
1358           double d2 = dist2( mpV, mp );
1359           dist2nID.insert( make_pair( d2, elem[j] ));
1360           avgDist2 += d2 / elem.GetNP();
1361           if ( !nbV )
1362             totNbSeg++, totSegLen+= sqrt( dist2( mp, ngMesh.Point( elem[(j+1)%elem.GetNP()])));
1363         }
1364         double dist = dist2nID.begin()->first; //avgDist2;
1365         if ( dist < closeDist2 )
1366           vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= dist;
1367       }
1368       dist2VData.insert( make_pair( closeDist2, vData ));
1369     }
1370
1371     if ( totNbSeg == 0 ) break;
1372     double avgSegLen = totSegLen / totNbSeg;
1373
1374     // Loop on vertices to add triangles
1375
1376     multimap< double, TIntVSoData >::iterator dist_vData = dist2VData.begin();
1377     for ( ; dist_vData != dist2VData.end(); ++dist_vData )
1378     {
1379       double closeDist2   = dist_vData->first;
1380       TIntVSoData & vData = dist_vData->second;
1381
1382       const netgen::MeshPoint& mpV = ngMesh.Point( vData.ngId );
1383
1384       // try to find more close face among ones added for internal vertices
1385       for (i = nbFaceInit+1; i <= ngMesh.GetNSE(); ++i)
1386       {
1387         double avgDist2 = 0;
1388         multimap< double, int> dist2nID;
1389         const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1390         for ( int j = 0; j < elem.GetNP(); ++j)
1391         {
1392           double d = dist2( mpV, ngMesh.Point( elem[j] ));
1393           dist2nID.insert( make_pair( d, elem[j] ));
1394           avgDist2 += d / elem.GetNP();
1395           if ( avgDist2 < closeDist2 )
1396             vData.ngIdClose= i, vData.ngIdCloseN= dist2nID.begin()->second, closeDist2= avgDist2;
1397         }
1398       }
1399       // sort nodes of the closest face by angle with vector from V to the closest node
1400       const double tol = numeric_limits<double>::min();
1401       map< double, int > angle2ID;
1402       const netgen::Element2d& closeFace = ngMesh.SurfaceElement( vData.ngIdClose );
1403       netgen::MeshPoint mp[2];
1404       mp[0] = ngMesh.Point( vData.ngIdCloseN );
1405       gp_XYZ p1( NGPOINT_COORDS( mp[0] ));
1406       gp_XYZ pV( NGPOINT_COORDS( mpV ));
1407       gp_Vec v2p1( pV, p1 );
1408       double distN1 = v2p1.Magnitude();
1409       if ( distN1 <= tol ) continue;
1410       v2p1 /= distN1;
1411       for ( int j = 0; j < closeFace.GetNP(); ++j)
1412       {
1413         mp[1] = ngMesh.Point( closeFace[j] );
1414         gp_Vec v2p( pV, gp_Pnt( NGPOINT_COORDS( mp[1] )) );
1415         angle2ID.insert( make_pair( v2p1.Angle( v2p ), closeFace[j]));
1416       }
1417       // get node with angle of 60 degrees or greater
1418       map< double, int >::iterator angle_id = angle2ID.lower_bound( 60. * M_PI / 180. );
1419       if ( angle_id == angle2ID.end() ) angle_id = --angle2ID.end();
1420       const double minAngle = 30. * M_PI / 180.;
1421       const double angle = angle_id->first;
1422       bool angleOK = ( angle > minAngle );
1423
1424       // find points to create a triangle
1425       netgen::Element2d tri(3);
1426       tri.SetIndex ( 1 );
1427       tri[0] = vData.ngId;
1428       tri[1] = vData.ngIdCloseN; // to use the closest nodes
1429       tri[2] = angle_id->second; // to use the node with best angle
1430
1431       // decide whether to use the closest node and the node with best angle or to create new ones
1432       for ( int isBestAngleN = 0; isBestAngleN < 2; ++isBestAngleN )
1433       {
1434         bool createNew = !angleOK, distOK = true;
1435         double distFromV;
1436         int triInd = isBestAngleN ? 2 : 1;
1437         mp[isBestAngleN] = ngMesh.Point( tri[triInd] );
1438         if ( isBestAngleN )
1439         {
1440           if ( angleOK )
1441           {
1442             double distN2 = sqrt( dist2( mpV, mp[isBestAngleN]));
1443             createNew = ( fabs( distN2 - distN1 ) > 0.25 * distN1 );
1444           }
1445           else if ( angle < tol )
1446           {
1447             v2p1.SetX( v2p1.X() + 1e-3 );
1448           }
1449           distFromV = distN1;
1450         }
1451         else
1452         {
1453           double segLenHint = ngMesh.GetH( ngMesh.Point( vData.ngId ));
1454           bool avgLenOK  = ( avgSegLen < 0.75 * distN1 );
1455           bool hintLenOK = ( segLenHint  < 0.75 * distN1 );
1456           createNew = (createNew || avgLenOK || hintLenOK );
1457           // we create a new node not closer than 0.5 to the closest face
1458           // in order not to clash with other close face
1459           double r = min( 0.5, ( hintLenOK ? segLenHint : avgSegLen ) / distN1 );
1460           distFromV = r * distN1;
1461         }
1462         if ( createNew )
1463         {
1464           // create a new point, between the node and the vertex if angleOK
1465           gp_XYZ p( NGPOINT_COORDS( mp[isBestAngleN] ));
1466           gp_Vec v2p( pV, p ); v2p.Normalize();
1467           if ( isBestAngleN && !angleOK )
1468             p = p1 + gp_Dir( v2p.XYZ() - v2p1.XYZ()).XYZ() * distN1 * 0.95;
1469           else
1470             p = pV + v2p.XYZ() * distFromV;
1471
1472           if ( !isBestAngleN ) p1 = p, distN1 = distFromV;
1473
1474           mp[isBestAngleN].SetPoint( netgen::Point<3> (p.X(), p.Y(), p.Z()));
1475           ngMesh.AddPoint ( mp[isBestAngleN], 1, netgen::SURFACEPOINT );
1476           tri[triInd] = ngMesh.GetNP();
1477           nodeVec.push_back( helper.AddNode( p.X(), p.Y(), p.Z()) );
1478         }
1479       }
1480       ngMesh.AddSurfaceElement (tri);
1481       swap( tri[1], tri[2] );
1482       ngMesh.AddSurfaceElement (tri);
1483
1484 #ifdef DUMP_TRIANGLES_SCRIPT
1485       py << "n1 = m.AddNode( "<< mpV.X()<<", "<< mpV.Y()<<", "<< mpV.Z()<<") "<< endl
1486          << "n2 = m.AddNode( "<< mp[0].X()<<", "<< mp[0].Y()<<", "<< mp[0].Z()<<") "<< endl
1487          << "n3 = m.AddNode( "<< mp[1].X()<<", "<< mp[1].Y()<<", "<< mp[1].Z()<<" )" << endl
1488          << "m.AddFace([n1,n2,n3])" << endl;
1489 #endif
1490     } // loop on internal vertices of a solid
1491
1492   } // loop on solids with internal vertices
1493 }
1494
1495 //================================================================================
1496 /*!
1497  * \brief Fill netgen mesh with segments of a FACE
1498  *  \param ngMesh - netgen mesh
1499  *  \param geom - container of OCCT geometry to mesh
1500  *  \param wires - data of nodes on FACE boundary
1501  *  \param helper - mesher helper holding the FACE
1502  *  \param nodeVec - vector of nodes in which node index == netgen ID
1503  *  \retval SMESH_ComputeErrorPtr - error description 
1504  */
1505 //================================================================================
1506
1507 SMESH_ComputeErrorPtr
1508 NETGENPlugin_Mesher::AddSegmentsToMesh(netgen::Mesh&                    ngMesh,
1509                                        netgen::OCCGeometry&             geom,
1510                                        const TSideVector&               wires,
1511                                        SMESH_MesherHelper&              helper,
1512                                        vector< const SMDS_MeshNode* > & nodeVec)
1513 {
1514   // ----------------------------
1515   // Check wires and count nodes
1516   // ----------------------------
1517   int nbNodes = 0;
1518   for ( int iW = 0; iW < wires.size(); ++iW )
1519   {
1520     StdMeshers_FaceSidePtr wire = wires[ iW ];
1521     if ( wire->MissVertexNode() )
1522     {
1523       // Commented for issue 0020960. It worked for the case, let's wait for case where it doesn't.
1524       // It seems that there is no reason for this limitation
1525 //       return TError
1526 //         (new SMESH_ComputeError(COMPERR_BAD_INPUT_MESH, "Missing nodes on vertices"));
1527     }
1528     const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1529     if ( uvPtVec.size() != wire->NbPoints() )
1530       return SMESH_ComputeError::New(COMPERR_BAD_INPUT_MESH,
1531                                      SMESH_Comment("Unexpected nb of points on wire ") << iW
1532                                      << ": " << uvPtVec.size()<<" != "<<wire->NbPoints());
1533     nbNodes += wire->NbPoints();
1534   }
1535   nodeVec.reserve( nodeVec.size() + nbNodes + 1 );
1536   if ( nodeVec.empty() )
1537     nodeVec.push_back( 0 );
1538
1539   // -----------------
1540   // Fill netgen mesh
1541   // -----------------
1542
1543   const bool wasNgMeshEmpty = ( ngMesh.GetNP() < 1 ); /* true => this method is called by
1544                                                          NETGENPlugin_NETGEN_2D_ONLY */
1545
1546   // map for nodes on vertices since they can be shared between wires
1547   // ( issue 0020676, face_int_box.brep) and nodes built by NETGEN
1548   map<const SMDS_MeshNode*, int > node2ngID;
1549   if ( !wasNgMeshEmpty ) // fill node2ngID with nodes built by NETGEN
1550   {
1551     set< int > subIDs; // ids of sub-shapes of the FACE
1552     for ( int iW = 0; iW < wires.size(); ++iW )
1553     {
1554       StdMeshers_FaceSidePtr wire = wires[ iW ];
1555       for ( int iE = 0, nbE = wire->NbEdges(); iE < nbE; ++iE )
1556       {
1557         subIDs.insert( wire->EdgeID( iE ));
1558         subIDs.insert( helper.GetMeshDS()->ShapeToIndex( wire->FirstVertex( iE )));
1559       }
1560     }
1561     for ( size_t ngID = 1; ngID < nodeVec.size(); ++ngID )
1562       if ( subIDs.count( nodeVec[ngID]->getshapeId() ))
1563         node2ngID.insert( make_pair( nodeVec[ngID], ngID ));
1564   }
1565
1566   const int solidID = 0, faceID = geom.fmap.FindIndex( helper.GetSubShape() );
1567   if ( ngMesh.GetNFD() < 1 )
1568     ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(faceID, solidID, solidID, 0));
1569
1570   for ( int iW = 0; iW < wires.size(); ++iW )
1571   {
1572     StdMeshers_FaceSidePtr wire = wires[ iW ];
1573     const vector<UVPtStruct>& uvPtVec = wire->GetUVPtStruct();
1574     const int nbSegments = wire->NbPoints() - 1;
1575
1576     // assure the 1st node to be in node2ngID, which is needed to correctly
1577     // "close chain of segments" (see below) in case if the 1st node is not
1578     // onVertex because it is on a Viscous layer
1579     node2ngID.insert( make_pair( uvPtVec[ 0 ].node, ngMesh.GetNP() + 1 ));
1580
1581     // compute length of every segment
1582     vector<double> segLen( nbSegments );
1583     for ( int i = 0; i < nbSegments; ++i )
1584       segLen[i] = SMESH_TNodeXYZ( uvPtVec[ i ].node ).Distance( uvPtVec[ i+1 ].node );
1585
1586     int edgeID = 1, posID = -2;
1587     bool isInternalWire = false;
1588     double vertexNormPar = 0;
1589     const int prevNbNGSeg = ngMesh.GetNSeg();
1590     for ( int i = 0; i < nbSegments; ++i ) // loop on segments
1591     {
1592       // Add the first point of a segment
1593
1594       const SMDS_MeshNode * n = uvPtVec[ i ].node;
1595       const int posShapeID = n->getshapeId();
1596       bool onVertex = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_VERTEX );
1597       bool onEdge   = ( n->GetPosition()->GetTypeOfPosition() == SMDS_TOP_EDGE   );
1598
1599       // skip nodes on degenerated edges
1600       if ( helper.IsDegenShape( posShapeID ) &&
1601            helper.IsDegenShape( uvPtVec[ i+1 ].node->getshapeId() ))
1602         continue;
1603
1604       int ngID1 = ngMesh.GetNP() + 1, ngID2 = ngID1+1;
1605       if ( onVertex || ( !wasNgMeshEmpty && onEdge ))
1606         ngID1 = node2ngID.insert( make_pair( n, ngID1 )).first->second;
1607       if ( ngID1 > ngMesh.GetNP() )
1608       {
1609         netgen::MeshPoint mp( netgen::Point<3> (n->X(), n->Y(), n->Z()) );
1610         ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1611         nodeVec.push_back( n );
1612       }
1613       else // n is in ngMesh already, and ngID2 in prev segment is wrong
1614       {
1615         ngID2 = ngMesh.GetNP() + 1;
1616         if ( i > 0 ) // prev segment belongs to same wire
1617         {
1618           netgen::Segment& prevSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
1619           prevSeg[1] = ngID1;
1620         }
1621       }
1622
1623       // Add the segment
1624
1625       netgen::Segment seg;
1626
1627       seg[0]     = ngID1;                // ng node id
1628       seg[1]     = ngID2;                // ng node id
1629       seg.edgenr = ngMesh.GetNSeg() + 1; // ng segment id
1630       seg.si     = faceID;               // = geom.fmap.FindIndex (face);
1631
1632       for ( int iEnd = 0; iEnd < 2; ++iEnd)
1633       {
1634         const UVPtStruct& pnt = uvPtVec[ i + iEnd ];
1635
1636         seg.epgeominfo[ iEnd ].dist = pnt.param; // param on curve
1637         seg.epgeominfo[ iEnd ].u    = pnt.u;
1638         seg.epgeominfo[ iEnd ].v    = pnt.v;
1639
1640         // find out edge id and node parameter on edge
1641         onVertex = ( pnt.normParam + 1e-10 > vertexNormPar );
1642         if ( onVertex || posShapeID != posID )
1643         {
1644           // get edge id
1645           double normParam = pnt.normParam;
1646           if ( onVertex )
1647             normParam = 0.5 * ( uvPtVec[ i ].normParam + uvPtVec[ i+1 ].normParam );
1648           int edgeIndexInWire = wire->EdgeIndex( normParam );
1649           vertexNormPar = wire->LastParameter( edgeIndexInWire );
1650           const TopoDS_Edge& edge = wire->Edge( edgeIndexInWire );
1651           edgeID = geom.emap.FindIndex( edge );
1652           posID  = posShapeID;
1653           isInternalWire = ( edge.Orientation() == TopAbs_INTERNAL );
1654           // if ( onVertex ) // param on curve is different on each of two edges
1655           //   seg.epgeominfo[ iEnd ].dist = helper.GetNodeU( edge, pnt.node );
1656         }
1657         seg.epgeominfo[ iEnd ].edgenr = edgeID; //  = geom.emap.FindIndex(edge);
1658       }
1659
1660       ngMesh.AddSegment (seg);
1661       {
1662         // restrict size of elements near the segment
1663         SMESH_TNodeXYZ np1( n ), np2( uvPtVec[ i+1 ].node );
1664         // get an average size of adjacent segments to avoid sharp change of
1665         // element size (regression on issue 0020452, note 0010898)
1666         int iPrev = SMESH_MesherHelper::WrapIndex( i-1, nbSegments );
1667         int iNext = SMESH_MesherHelper::WrapIndex( i+1, nbSegments );
1668         double avgH = ( segLen[ iPrev ] + segLen[ i ] + segLen[ iNext ]) / 3;
1669
1670         RestrictLocalSize( ngMesh, 0.5*(np1+np2), avgH );
1671       }
1672       if ( isInternalWire )
1673       {
1674         swap (seg[0], seg[1]);
1675         swap( seg.epgeominfo[0], seg.epgeominfo[1] );
1676         seg.edgenr = ngMesh.GetNSeg() + 1; // segment id
1677         ngMesh.AddSegment (seg);
1678       }
1679     } // loop on segments on a wire
1680
1681     // close chain of segments
1682     if ( nbSegments > 0 )
1683     {
1684       netgen::Segment& lastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() - int( isInternalWire));
1685       const SMDS_MeshNode * lastNode = uvPtVec.back().node;
1686       lastSeg[1] = node2ngID.insert( make_pair( lastNode, lastSeg[1] )).first->second;
1687       if ( lastSeg[1] > ngMesh.GetNP() )
1688       {
1689         netgen::MeshPoint mp( netgen::Point<3> (lastNode->X(), lastNode->Y(), lastNode->Z()) );
1690         ngMesh.AddPoint ( mp, 1, netgen::EDGEPOINT );
1691         nodeVec.push_back( lastNode );
1692       }
1693       if ( isInternalWire )
1694       {
1695         netgen::Segment& realLastSeg = ngMesh.LineSegment( ngMesh.GetNSeg() );
1696         realLastSeg[0] = lastSeg[1];
1697       }
1698     }
1699
1700 #ifdef DUMP_SEGMENTS
1701     cout << "BEGIN WIRE " << iW << endl;
1702     for ( int i = prevNbNGSeg+1; i <= ngMesh.GetNSeg(); ++i )
1703     {
1704       netgen::Segment& seg = ngMesh.LineSegment( i );
1705       if ( i > 1 ) {
1706         netgen::Segment& prevSeg = ngMesh.LineSegment( i-1 );
1707         if ( seg[0] == prevSeg[1] && seg[1] == prevSeg[0] )
1708         {
1709           cout << "Segment: " << seg.edgenr << endl << "\tis REVRESE of the previous one" << endl;
1710           continue;
1711         }
1712       }
1713       cout << "Segment: " << seg.edgenr << endl
1714            << "\tp1: " << seg[0] << endl
1715            << "\tp2: " << seg[1] << endl
1716            << "\tp0 param: " << seg.epgeominfo[ 0 ].dist << endl
1717            << "\tp0 uv: " << seg.epgeominfo[ 0 ].u <<", "<< seg.epgeominfo[ 0 ].v << endl
1718            << "\tp0 edge: " << seg.epgeominfo[ 0 ].edgenr << endl
1719            << "\tp1 param: " << seg.epgeominfo[ 1 ].dist << endl
1720            << "\tp1 uv: " << seg.epgeominfo[ 1 ].u <<", "<< seg.epgeominfo[ 1 ].v << endl
1721            << "\tp1 edge: " << seg.epgeominfo[ 1 ].edgenr << endl;
1722     }
1723     cout << "--END WIRE " << iW << endl;
1724 #endif
1725
1726   } // loop on WIREs of a FACE
1727
1728   // add a segment instead of an internal vertex
1729   if ( wasNgMeshEmpty )
1730   {
1731     NETGENPlugin_Internals intShapes( *helper.GetMesh(), helper.GetSubShape(), /*is3D=*/false );
1732     AddIntVerticesInFaces( geom, ngMesh, nodeVec, intShapes );
1733   }
1734   ngMesh.CalcSurfacesOfNode();
1735
1736   return TError();
1737 }
1738
1739 //================================================================================
1740 /*!
1741  * \brief Fill SMESH mesh according to contents of netgen mesh
1742  *  \param occgeo - container of OCCT geometry to mesh
1743  *  \param ngMesh - netgen mesh
1744  *  \param initState - bn of entities in netgen mesh before computing
1745  *  \param sMesh - SMESH mesh to fill in
1746  *  \param nodeVec - vector of nodes in which node index == netgen ID
1747  *  \retval int - error
1748  */
1749 //================================================================================
1750
1751 int NETGENPlugin_Mesher::FillSMesh(const netgen::OCCGeometry&          occgeo,
1752                                    netgen::Mesh&                       ngMesh,
1753                                    const NETGENPlugin_ngMeshInfo&      initState,
1754                                    SMESH_Mesh&                         sMesh,
1755                                    std::vector<const SMDS_MeshNode*>&  nodeVec,
1756                                    SMESH_Comment&                      comment)
1757 {
1758   int nbNod = ngMesh.GetNP();
1759   int nbSeg = ngMesh.GetNSeg();
1760   int nbFac = ngMesh.GetNSE();
1761   int nbVol = ngMesh.GetNE();
1762
1763   SMESHDS_Mesh* meshDS = sMesh.GetMeshDS();
1764
1765   // -------------------------------------
1766   // Create and insert nodes into nodeVec
1767   // -------------------------------------
1768
1769   nodeVec.resize( nbNod + 1 );
1770   int i, nbInitNod = initState._nbNodes;
1771   for (i = nbInitNod+1; i <= nbNod; ++i )
1772   {
1773     const netgen::MeshPoint& ngPoint = ngMesh.Point(i);
1774     SMDS_MeshNode* node = NULL;
1775     TopoDS_Vertex aVert;
1776     // First, netgen creates nodes on vertices in occgeo.vmap,
1777     // so node index corresponds to vertex index
1778     // but (issue 0020776) netgen does not create nodes with equal coordinates
1779     if ( i-nbInitNod <= occgeo.vmap.Extent() )
1780     {
1781       gp_Pnt p ( NGPOINT_COORDS(ngPoint) );
1782       for (int iV = i-nbInitNod; aVert.IsNull() && iV <= occgeo.vmap.Extent(); ++iV)
1783       {
1784         aVert = TopoDS::Vertex( occgeo.vmap( iV ) );
1785         gp_Pnt pV = BRep_Tool::Pnt( aVert );
1786         if ( p.SquareDistance( pV ) > 1e-20 )
1787           aVert.Nullify();
1788         else
1789           node = const_cast<SMDS_MeshNode*>( SMESH_Algo::VertexNode( aVert, meshDS ));
1790       }
1791     }
1792     if (!node) // node not found on vertex
1793     {
1794       node = meshDS->AddNode( NGPOINT_COORDS( ngPoint ));
1795       if (!aVert.IsNull())
1796         meshDS->SetNodeOnVertex(node, aVert);
1797     }
1798     nodeVec[i] = node;
1799   }
1800
1801   // -------------------------------------------
1802   // Create mesh segments along geometric edges
1803   // -------------------------------------------
1804
1805   int nbInitSeg = initState._nbSegments;
1806   for (i = nbInitSeg+1; i <= nbSeg; ++i )
1807   {
1808     const netgen::Segment& seg = ngMesh.LineSegment(i);
1809     TopoDS_Edge aEdge;
1810     int pinds[3] = { seg.pnums[0], seg.pnums[1], seg.pnums[2] };
1811     int nbp = 0;
1812     double param2 = 0;
1813     for (int j=0; j < 3; ++j)
1814     {
1815       int pind = pinds[j];
1816       if (pind <= 0 || !nodeVec_ACCESS(pind))
1817         break;
1818       ++nbp;
1819       double param;
1820       if (j < 2)
1821       {
1822         if (aEdge.IsNull())
1823         {
1824           int aGeomEdgeInd = seg.epgeominfo[j].edgenr;
1825           if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
1826             aEdge = TopoDS::Edge(occgeo.emap(aGeomEdgeInd));
1827         }
1828         param = seg.epgeominfo[j].dist;
1829         param2 += param;
1830       }
1831       else // middle point
1832       {
1833         param = param2 * 0.5;
1834       }
1835       if (!aEdge.IsNull() && nodeVec_ACCESS(pind)->getshapeId() < 1)
1836       {
1837         meshDS->SetNodeOnEdge(nodeVec_ACCESS(pind), aEdge, param);
1838       }
1839     }
1840     if ( nbp > 1 )
1841     {
1842       SMDS_MeshEdge* edge = 0;
1843       if (nbp == 2) // second order ?
1844       {
1845         if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1])))
1846           continue;
1847         edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]));
1848       }
1849       else
1850       {
1851         if ( meshDS->FindEdge( nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
1852                                nodeVec_ACCESS(pinds[2])))
1853           continue;
1854         edge = meshDS->AddEdge(nodeVec_ACCESS(pinds[0]), nodeVec_ACCESS(pinds[1]),
1855                                nodeVec_ACCESS(pinds[2]));
1856       }
1857       if (!edge)
1858       {
1859         if ( comment.empty() ) comment << "Cannot create a mesh edge";
1860         MESSAGE("Cannot create a mesh edge");
1861         nbSeg = nbFac = nbVol = 0;
1862         break;
1863       }
1864       if ( !aEdge.IsNull() && edge->getshapeId() < 1 )
1865         meshDS->SetMeshElementOnShape(edge, aEdge);
1866     }
1867     else if ( comment.empty() )
1868     {
1869       comment << "Invalid netgen segment #" << i;
1870     }
1871   }
1872
1873   // ----------------------------------------
1874   // Create mesh faces along geometric faces
1875   // ----------------------------------------
1876
1877   int nbInitFac = initState._nbFaces;
1878   int quadFaceID = ngMesh.GetNFD() + 1;
1879   if ( nbInitFac < nbFac )
1880     // add a faces descriptor to exclude qudrangle elements generated by NETGEN
1881     // from computation of 3D mesh
1882     ngMesh.AddFaceDescriptor (netgen::FaceDescriptor(quadFaceID, /*solid1=*/0, /*solid2=*/0, 0));
1883
1884   for (i = nbInitFac+1; i <= nbFac; ++i )
1885   {
1886     const netgen::Element2d& elem = ngMesh.SurfaceElement(i);
1887     int aGeomFaceInd = elem.GetIndex();
1888     TopoDS_Face aFace;
1889     if (aGeomFaceInd > 0 && aGeomFaceInd <= occgeo.fmap.Extent())
1890       aFace = TopoDS::Face(occgeo.fmap(aGeomFaceInd));
1891     vector<SMDS_MeshNode*> nodes;
1892     for (int j=1; j <= elem.GetNP(); ++j)
1893     {
1894       int pind = elem.PNum(j);
1895       if ( pind < 1 || pind >= nodeVec.size() )
1896         break;
1897       if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind))
1898       {
1899         nodes.push_back(node);
1900         if (!aFace.IsNull() && node->getshapeId() < 1)
1901         {
1902           const netgen::PointGeomInfo& pgi = elem.GeomInfoPi(j);
1903           meshDS->SetNodeOnFace(node, aFace, pgi.u, pgi.v);
1904         }
1905       }
1906     }
1907     if ( nodes.size() != elem.GetNP() )
1908     {
1909       if ( comment.empty() )
1910         comment << "Invalid netgen 2d element #" << i;
1911       continue; // bad node ids
1912     }
1913     SMDS_MeshFace* face = NULL;
1914     switch (elem.GetType())
1915     {
1916     case netgen::TRIG:
1917       face = meshDS->AddFace(nodes[0],nodes[1],nodes[2]);
1918       break;
1919     case netgen::QUAD:
1920       face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3]);
1921       // exclude qudrangle elements from computation of 3D mesh
1922       const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
1923       break;
1924     case netgen::TRIG6:
1925       face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[5],nodes[3],nodes[4]);
1926       break;
1927     case netgen::QUAD8:
1928       face = meshDS->AddFace(nodes[0],nodes[1],nodes[2],nodes[3],
1929                              nodes[4],nodes[7],nodes[5],nodes[6]);
1930       // exclude qudrangle elements from computation of 3D mesh
1931       const_cast< netgen::Element2d& >( elem ).SetIndex( quadFaceID );
1932       break;
1933     default:
1934       MESSAGE("NETGEN created a face of unexpected type, ignoring");
1935       continue;
1936     }
1937     if (!face)
1938     {
1939       if ( comment.empty() ) comment << "Cannot create a mesh face";
1940       MESSAGE("Cannot create a mesh face");
1941       nbSeg = nbFac = nbVol = 0;
1942       break;
1943     }
1944     if (!aFace.IsNull())
1945       meshDS->SetMeshElementOnShape(face, aFace);
1946   }
1947
1948   // ------------------
1949   // Create tetrahedra
1950   // ------------------
1951
1952   for (i = 1; i <= nbVol; ++i)
1953   {
1954     const netgen::Element& elem = ngMesh.VolumeElement(i);      
1955     int aSolidInd = elem.GetIndex();
1956     TopoDS_Solid aSolid;
1957     if (aSolidInd > 0 && aSolidInd <= occgeo.somap.Extent())
1958       aSolid = TopoDS::Solid(occgeo.somap(aSolidInd));
1959     vector<SMDS_MeshNode*> nodes;
1960     for (int j=1; j <= elem.GetNP(); ++j)
1961     {
1962       int pind = elem.PNum(j);
1963       if ( pind < 1 || pind >= nodeVec.size() )
1964         break;
1965       if ( SMDS_MeshNode* node = nodeVec_ACCESS(pind) )
1966       {
1967         nodes.push_back(node);
1968         if ( !aSolid.IsNull() && node->getshapeId() < 1 )
1969           meshDS->SetNodeInVolume(node, aSolid);
1970       }
1971     }
1972     if ( nodes.size() != elem.GetNP() )
1973     {
1974       if ( comment.empty() )
1975         comment << "Invalid netgen 3d element #" << i;
1976       continue;
1977     }
1978     SMDS_MeshVolume* vol = NULL;
1979     switch (elem.GetType())
1980     {
1981     case netgen::TET:
1982       vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3]);
1983       break;
1984     case netgen::TET10:
1985       vol = meshDS->AddVolume(nodes[0],nodes[1],nodes[2],nodes[3],
1986                               nodes[4],nodes[7],nodes[5],nodes[6],nodes[8],nodes[9]);
1987       break;
1988     default:
1989       MESSAGE("NETGEN created a volume of unexpected type, ignoring");
1990       continue;
1991     }
1992     if (!vol)
1993     {
1994       if ( comment.empty() ) comment << "Cannot create a mesh volume";
1995       MESSAGE("Cannot create a mesh volume");
1996       nbSeg = nbFac = nbVol = 0;
1997       break;
1998     }
1999     if (!aSolid.IsNull())
2000       meshDS->SetMeshElementOnShape(vol, aSolid);
2001   }
2002   return comment.empty() ? 0 : 1;
2003 }
2004
2005 namespace
2006 {
2007   //================================================================================
2008   /*!
2009    * \brief Restrict size of elements on the given edge 
2010    */
2011   //================================================================================
2012
2013   void setLocalSize(const TopoDS_Edge& edge,
2014                     double             size,
2015                     netgen::Mesh&      mesh)
2016   {
2017     const int nb = 1000;
2018     Standard_Real u1, u2;
2019     Handle(Geom_Curve) curve = BRep_Tool::Curve(edge, u1, u2);
2020     if ( curve.IsNull() )
2021     {
2022       TopoDS_Iterator vIt( edge );
2023       if ( !vIt.More() ) return;
2024       gp_Pnt p = BRep_Tool::Pnt( TopoDS::Vertex( vIt.Value() ));
2025       NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
2026     }
2027     else
2028     {
2029       Standard_Real delta = (u2-u1)/nb;
2030       for(int i=0; i<nb; i++)
2031       {
2032         Standard_Real u = u1 + delta*i;
2033         gp_Pnt p = curve->Value(u);
2034         NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), size );
2035         netgen::Point3d pi(p.X(), p.Y(), p.Z());
2036         double resultSize = mesh.GetH(pi);
2037         if ( resultSize - size > 0.1*size )
2038           // netgen does restriction iff oldH/newH > 1.2 (localh.cpp:136)
2039           NETGENPlugin_Mesher::RestrictLocalSize( mesh, p.XYZ(), resultSize/1.201 );
2040       }
2041     }
2042   }
2043
2044   //================================================================================
2045   /*!
2046    * \brief Convert error into text
2047    */
2048   //================================================================================
2049
2050   std::string text(int err)
2051   {
2052     if ( !err )
2053       return string("");
2054     return
2055       SMESH_Comment("Error in netgen::OCCGenerateMesh() at ") << netgen::multithread.task;
2056   }
2057
2058   //================================================================================
2059   /*!
2060    * \brief Convert exception into text
2061    */
2062   //================================================================================
2063
2064   std::string text(Standard_Failure& ex)
2065   {
2066     SMESH_Comment str("Exception in netgen::OCCGenerateMesh()");
2067     str << " at " << netgen::multithread.task
2068         << ": " << ex.DynamicType()->Name();
2069     if ( ex.GetMessageString() && strlen( ex.GetMessageString() ))
2070       str << ": " << ex.GetMessageString();
2071     return str;
2072   }
2073   //================================================================================
2074   /*!
2075    * \brief Convert exception into text
2076    */
2077   //================================================================================
2078
2079   std::string text(netgen::NgException& ex)
2080   {
2081     SMESH_Comment str("NgException");
2082     if ( strlen( netgen::multithread.task ) > 0 )
2083       str << " at " << netgen::multithread.task;
2084     str << ": " << ex.What();
2085     return str;
2086   }
2087
2088   const double edgeMeshingTime = 0.001;
2089   const double faceMeshingTime = 0.019;
2090   const double edgeFaceMeshingTime = edgeMeshingTime + faceMeshingTime;
2091   const double faceOptimizTime = 0.06;
2092   const double voluMeshingTime = 0.15;
2093   const double volOptimizeTime = 0.77;
2094 }
2095
2096 //=============================================================================
2097 /*!
2098  * Here we are going to use the NETGEN mesher
2099  */
2100 //=============================================================================
2101
2102 bool NETGENPlugin_Mesher::Compute()
2103 {
2104   NETGENPlugin_NetgenLibWrapper ngLib;
2105
2106   netgen::MeshingParameters& mparams = netgen::mparam;
2107   MESSAGE("Compute with:\n"
2108           " max size = " << mparams.maxh << "\n"
2109           " segments per edge = " << mparams.segmentsperedge);
2110   MESSAGE("\n"
2111           " growth rate = " << mparams.grading << "\n"
2112           " elements per radius = " << mparams.curvaturesafety << "\n"
2113           " second order = " << mparams.secondorder << "\n"
2114           " quad allowed = " << mparams.quad);
2115   //cout << " quad allowed = " << mparams.quad<<endl;
2116
2117   SMESH_ComputeErrorPtr error = SMESH_ComputeError::New();
2118
2119   static string debugFile = "/tmp/ngMesh.py"; /* to call toPython( ngMesh, debugFile )
2120                                                  while debugging netgen */
2121   // -------------------------
2122   // Prepare OCC geometry
2123   // -------------------------
2124
2125   netgen::OCCGeometry occgeo;
2126   list< SMESH_subMesh* > meshedSM[3]; // for 0-2 dimensions
2127   NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
2128   PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
2129   _occgeom = &occgeo;
2130
2131   _totalTime = edgeFaceMeshingTime;
2132   if ( _optimize )
2133     _totalTime += faceOptimizTime;
2134   if ( _isVolume )
2135     _totalTime += voluMeshingTime + ( _optimize ? volOptimizeTime : 0 );
2136   double doneTime = 0;
2137   _ticTime = -1;
2138   _progressTic = 1;
2139   _curShapeIndex = -1;
2140
2141   // -------------------------
2142   // Generate the mesh
2143   // -------------------------
2144
2145   _ngMesh = NULL;
2146   NETGENPlugin_ngMeshInfo initState; // it remembers size of ng mesh equal to size of Smesh
2147
2148   SMESH_Comment comment;
2149   int err = 0;
2150
2151   // vector of nodes in which node index == netgen ID
2152   vector< const SMDS_MeshNode* > nodeVec;
2153   
2154   {
2155     // ----------------
2156     // compute 1D mesh
2157     // ----------------
2158     if ( _simpleHyp )
2159     {
2160       // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
2161       mparams.uselocalh = false;
2162       mparams.grading = 0.8; // not limitited size growth
2163
2164       if ( _simpleHyp->GetNumberOfSegments() )
2165         // nb of segments
2166         mparams.maxh = occgeo.boundingbox.Diam();
2167       else
2168         // segment length
2169         mparams.maxh = _simpleHyp->GetLocalLength();
2170     }
2171
2172     if ( mparams.maxh == 0.0 )
2173       mparams.maxh = occgeo.boundingbox.Diam();
2174     if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
2175       mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
2176
2177     // Local size on faces
2178     occgeo.face_maxh = mparams.maxh;
2179
2180     // Let netgen create _ngMesh and calculate element size on not meshed shapes
2181 #ifndef NETGEN_V5
2182     char *optstr = 0;
2183 #endif
2184     int startWith = netgen::MESHCONST_ANALYSE;
2185     int endWith   = netgen::MESHCONST_ANALYSE;
2186     try
2187     {
2188       OCC_CATCH_SIGNALS;
2189 #ifdef NETGEN_V5
2190       err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2191 #else
2192       err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2193 #endif
2194       if(netgen::multithread.terminate)
2195         return false;
2196
2197       comment << text(err);
2198     }
2199     catch (Standard_Failure& ex)
2200     {
2201       comment << text(ex);
2202     }
2203     err = 0; //- MESHCONST_ANALYSE isn't so important step
2204     if ( !_ngMesh )
2205       return false;
2206     ngLib.setMesh(( Ng_Mesh*) _ngMesh );
2207
2208     _ngMesh->ClearFaceDescriptors(); // we make descriptors our-self
2209
2210     if ( _simpleHyp )
2211     {
2212       // Pass 1D simple parameters to NETGEN
2213       // --------------------------------
2214       int      nbSeg = _simpleHyp->GetNumberOfSegments();
2215       double segSize = _simpleHyp->GetLocalLength();
2216       for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
2217       {
2218         const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
2219         if ( nbSeg )
2220           segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
2221         setLocalSize( e, segSize, *_ngMesh );
2222       }
2223     }
2224     else // if ( ! _simpleHyp )
2225     {
2226       // Local size on vertices and edges
2227       // --------------------------------
2228       for(std::map<int,double>::const_iterator it=EdgeId2LocalSize.begin(); it!=EdgeId2LocalSize.end(); it++)
2229       {
2230         int key = (*it).first;
2231         double hi = (*it).second;
2232         const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2233         const TopoDS_Edge& e = TopoDS::Edge(shape);
2234         setLocalSize( e, hi, *_ngMesh );
2235       }
2236       for(std::map<int,double>::const_iterator it=VertexId2LocalSize.begin(); it!=VertexId2LocalSize.end(); it++)
2237       {
2238         int key = (*it).first;
2239         double hi = (*it).second;
2240         const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2241         const TopoDS_Vertex& v = TopoDS::Vertex(shape);
2242         gp_Pnt p = BRep_Tool::Pnt(v);
2243         NETGENPlugin_Mesher::RestrictLocalSize( *_ngMesh, p.XYZ(), hi );
2244       }
2245       for(map<int,double>::const_iterator it=FaceId2LocalSize.begin();
2246           it!=FaceId2LocalSize.end(); it++)
2247       {
2248         int key = (*it).first;
2249         double val = (*it).second;
2250         const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2251         int faceNgID = occgeo.fmap.FindIndex(shape);
2252         occgeo.SetFaceMaxH(faceNgID, val);
2253         for ( TopExp_Explorer edgeExp( shape, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
2254           setLocalSize( TopoDS::Edge( edgeExp.Current() ), val, *_ngMesh );
2255       }
2256     }
2257
2258     // Precompute internal edges (issue 0020676) in order to
2259     // add mesh on them correctly (twice) to netgen mesh
2260     if ( !err && internals.hasInternalEdges() )
2261     {
2262       // load internal shapes into OCCGeometry
2263       netgen::OCCGeometry intOccgeo;
2264       internals.getInternalEdges( intOccgeo.fmap, intOccgeo.emap, intOccgeo.vmap, meshedSM );
2265       intOccgeo.boundingbox = occgeo.boundingbox;
2266       intOccgeo.shape = occgeo.shape;
2267       intOccgeo.face_maxh.SetSize(intOccgeo.fmap.Extent());
2268       intOccgeo.face_maxh = netgen::mparam.maxh;
2269       netgen::Mesh *tmpNgMesh = NULL;
2270       try
2271       {
2272         OCC_CATCH_SIGNALS;
2273         // compute local H on internal shapes in the main mesh
2274         //OCCSetLocalMeshSize(intOccgeo, *_ngMesh); it deletes _ngMesh->localH
2275
2276         // let netgen create a temporary mesh
2277 #ifdef NETGEN_V5
2278         netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2279 #else
2280         netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2281 #endif
2282         if(netgen::multithread.terminate)
2283           return false;
2284
2285         // copy LocalH from the main to temporary mesh
2286         initState.transferLocalH( _ngMesh, tmpNgMesh );
2287
2288         // compute mesh on internal edges
2289         startWith = endWith = netgen::MESHCONST_MESHEDGES;
2290 #ifdef NETGEN_V5
2291         err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, mparams, startWith, endWith);
2292 #else
2293         err = netgen::OCCGenerateMesh(intOccgeo, tmpNgMesh, startWith, endWith, optstr);
2294 #endif
2295         comment << text(err);
2296       }
2297       catch (Standard_Failure& ex)
2298       {
2299         comment << text(ex);
2300         err = 1;
2301       }
2302       initState.restoreLocalH( tmpNgMesh );
2303
2304       // fill SMESH by netgen mesh
2305       vector< const SMDS_MeshNode* > tmpNodeVec;
2306       FillSMesh( intOccgeo, *tmpNgMesh, initState, *_mesh, tmpNodeVec, comment );
2307       err = ( err || !comment.empty() );
2308
2309       nglib::Ng_DeleteMesh((nglib::Ng_Mesh*)tmpNgMesh);
2310     }
2311
2312     // Fill _ngMesh with nodes and segments of computed submeshes
2313     if ( !err )
2314     {
2315       err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_0D ]) &&
2316                 FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_1D ]));
2317     }
2318     initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2319
2320     // Compute 1d mesh
2321     if (!err)
2322     {
2323       startWith = endWith = netgen::MESHCONST_MESHEDGES;
2324       try
2325       {
2326         OCC_CATCH_SIGNALS;
2327 #ifdef NETGEN_V5
2328         err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2329 #else
2330         err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2331 #endif
2332         if(netgen::multithread.terminate)
2333           return false;
2334
2335         comment << text(err);
2336       }
2337       catch (Standard_Failure& ex)
2338       {
2339         comment << text(ex);
2340         err = 1;
2341       }
2342     }
2343     if ( _isVolume )
2344       _ticTime = ( doneTime += edgeMeshingTime ) / _totalTime / _progressTic;
2345
2346     mparams.uselocalh = true; // restore as it is used at surface optimization
2347
2348     // ---------------------
2349     // compute surface mesh
2350     // ---------------------
2351     if (!err)
2352     {
2353       // Pass 2D simple parameters to NETGEN
2354       if ( _simpleHyp ) {
2355         if ( double area = _simpleHyp->GetMaxElementArea() ) {
2356           // face area
2357           mparams.maxh = sqrt(2. * area/sqrt(3.0));
2358           mparams.grading = 0.4; // moderate size growth
2359         }
2360         else {
2361           // length from edges
2362           if ( _ngMesh->GetNSeg() ) {
2363             double edgeLength = 0;
2364             TopTools_MapOfShape visitedEdges;
2365             for ( TopExp_Explorer exp( _shape, TopAbs_EDGE ); exp.More(); exp.Next() )
2366               if( visitedEdges.Add(exp.Current()) )
2367                 edgeLength += SMESH_Algo::EdgeLength( TopoDS::Edge( exp.Current() ));
2368             // we have to multiply length by 2 since for each TopoDS_Edge there
2369             // are double set of NETGEN edges, in other words, we have to
2370             // divide _ngMesh->GetNSeg() by 2.
2371             mparams.maxh = 2*edgeLength / _ngMesh->GetNSeg();
2372           }
2373           else {
2374             mparams.maxh = 1000;
2375           }
2376           mparams.grading = 0.2; // slow size growth
2377         }
2378         mparams.quad = _simpleHyp->GetAllowQuadrangles();
2379         mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2380         _ngMesh->SetGlobalH (mparams.maxh);
2381         netgen::Box<3> bb = occgeo.GetBoundingBox();
2382         bb.Increase (bb.Diam()/20);
2383         _ngMesh->SetLocalH (bb.PMin(), bb.PMax(), mparams.grading);
2384       }
2385
2386       // Care of vertices internal in faces (issue 0020676)
2387       if ( internals.hasInternalVertexInFace() )
2388       {
2389         // store computed segments in SMESH in order not to create SMESH
2390         // edges for ng segments added by AddIntVerticesInFaces()
2391         FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2392         // add segments to faces with internal vertices
2393         AddIntVerticesInFaces( occgeo, *_ngMesh, nodeVec, internals );
2394         initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2395       }
2396
2397       // Build viscous layers
2398       if ( _isViscousLayers2D )
2399       {
2400         if ( !internals.hasInternalVertexInFace() ) {
2401           FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2402           initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2403         }
2404         SMESH_ProxyMesh::Ptr viscousMesh;
2405         SMESH_MesherHelper   helper( *_mesh );
2406         for ( int faceID = 1; faceID <= occgeo.fmap.Extent(); ++faceID )
2407         {
2408           const TopoDS_Face& F = TopoDS::Face( occgeo.fmap( faceID ));
2409           viscousMesh = StdMeshers_ViscousLayers2D::Compute( *_mesh, F );
2410           if ( !viscousMesh )
2411             return false;
2412           // exclude from computation ng segments built on EDGEs of F
2413           for (int i = 1; i <= _ngMesh->GetNSeg(); i++)
2414           {
2415             netgen::Segment & seg = _ngMesh->LineSegment(i);
2416             if (seg.si == faceID)
2417               seg.si = 0;
2418           }
2419           // add new segments to _ngMesh instead of excluded ones
2420           helper.SetSubShape( F );
2421           TSideVector wires =
2422             StdMeshers_FaceSide::GetFaceWires( F, *_mesh, /*skipMediumNodes=*/true,
2423                                                error, viscousMesh );
2424           error = AddSegmentsToMesh( *_ngMesh, occgeo, wires, helper, nodeVec );
2425
2426           if ( !error ) error = SMESH_ComputeError::New();
2427         }
2428         initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2429       }
2430
2431       // Let netgen compute 2D mesh
2432       startWith = netgen::MESHCONST_MESHSURFACE;
2433       endWith = _optimize ? netgen::MESHCONST_OPTSURFACE : netgen::MESHCONST_MESHSURFACE;
2434       try
2435       {
2436         OCC_CATCH_SIGNALS;
2437 #ifdef NETGEN_V5
2438         err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2439 #else
2440         err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2441 #endif
2442         if(netgen::multithread.terminate)
2443           return false;
2444
2445         comment << text (err);
2446       }
2447       catch (Standard_Failure& ex)
2448       {
2449         comment << text(ex);
2450         //err = 1; -- try to make volumes anyway
2451       }
2452       catch (netgen::NgException exc)
2453       {
2454         comment << text(exc);
2455         //err = 1; -- try to make volumes anyway
2456       }
2457     }
2458     if ( _isVolume )
2459     {
2460       doneTime += faceMeshingTime + ( _optimize ? faceOptimizTime : 0 );
2461       _ticTime = doneTime / _totalTime / _progressTic;
2462     }
2463     // ---------------------
2464     // generate volume mesh
2465     // ---------------------
2466     // Fill _ngMesh with nodes and faces of computed 2D submeshes
2467     if ( !err && _isVolume && ( !meshedSM[ MeshDim_2D ].empty() || mparams.quad ))
2468     {
2469       // load SMESH with computed segments and faces
2470       FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2471
2472       // compute pyramids on quadrangles
2473       SMESH_ProxyMesh::Ptr proxyMesh;
2474       if ( _mesh->NbQuadrangles() > 0 )
2475         for ( int iS = 1; iS <= occgeo.somap.Extent(); ++iS )
2476         {
2477           StdMeshers_QuadToTriaAdaptor* Adaptor = new StdMeshers_QuadToTriaAdaptor;
2478           proxyMesh.reset( Adaptor );
2479
2480           int nbPyrams = _mesh->NbPyramids();
2481           Adaptor->Compute( *_mesh, occgeo.somap(iS) );
2482           if ( nbPyrams != _mesh->NbPyramids() )
2483           {
2484             list< SMESH_subMesh* > quadFaceSM;
2485             for (TopExp_Explorer face(occgeo.somap(iS), TopAbs_FACE); face.More(); face.Next())
2486               if ( Adaptor->GetProxySubMesh( face.Current() ))
2487               {
2488                 quadFaceSM.push_back( _mesh->GetSubMesh( face.Current() ));
2489                 meshedSM[ MeshDim_2D ].remove( quadFaceSM.back() );
2490               }
2491             FillNgMesh(occgeo, *_ngMesh, nodeVec, quadFaceSM, proxyMesh);
2492           }
2493         }
2494       // fill _ngMesh with faces of sub-meshes
2495       err = ! ( FillNgMesh(occgeo, *_ngMesh, nodeVec, meshedSM[ MeshDim_2D ]));
2496       initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2497       //toPython( _ngMesh, "/tmp/ngPython.py");
2498     }
2499     if (!err && _isVolume)
2500     {
2501       // Pass 3D simple parameters to NETGEN
2502       const NETGENPlugin_SimpleHypothesis_3D* simple3d =
2503         dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
2504       if ( simple3d ) {
2505         if ( double vol = simple3d->GetMaxElementVolume() ) {
2506           // max volume
2507           mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
2508           mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2509         }
2510         else {
2511           // length from faces
2512           mparams.maxh = _ngMesh->AverageH();
2513         }
2514         _ngMesh->SetGlobalH (mparams.maxh);
2515         mparams.grading = 0.4;
2516 #ifdef NETGEN_V5
2517         _ngMesh->CalcLocalH(mparams.grading);
2518 #else
2519         _ngMesh->CalcLocalH();
2520 #endif
2521       }
2522       // Care of vertices internal in solids and internal faces (issue 0020676)
2523       if ( internals.hasInternalVertexInSolid() || internals.hasInternalFaces() )
2524       {
2525         // store computed faces in SMESH in order not to create SMESH
2526         // faces for ng faces added here
2527         FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment );
2528         // add ng faces to solids with internal vertices
2529         AddIntVerticesInSolids( occgeo, *_ngMesh, nodeVec, internals );
2530         // duplicate mesh faces on internal faces
2531         FixIntFaces( occgeo, *_ngMesh, internals );
2532         initState = NETGENPlugin_ngMeshInfo(_ngMesh);
2533       }
2534       // Let netgen compute 3D mesh
2535       startWith = endWith = netgen::MESHCONST_MESHVOLUME;
2536       try
2537       {
2538         OCC_CATCH_SIGNALS;
2539 #ifdef NETGEN_V5
2540         err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2541 #else
2542         err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2543 #endif
2544         if(netgen::multithread.terminate)
2545           return false;
2546
2547         if ( comment.empty() ) // do not overwrite a previos error
2548           comment << text(err);
2549       }
2550       catch (Standard_Failure& ex)
2551       {
2552         if ( comment.empty() ) // do not overwrite a previos error
2553           comment << text(ex);
2554         err = 1;
2555       }
2556       catch (netgen::NgException exc)
2557       {
2558         if ( comment.empty() ) // do not overwrite a previos error
2559           comment << text(exc);
2560         err = 1;
2561       }
2562       _ticTime = ( doneTime += voluMeshingTime ) / _totalTime / _progressTic;
2563
2564       // Let netgen optimize 3D mesh
2565       if ( !err && _optimize )
2566       {
2567         startWith = endWith = netgen::MESHCONST_OPTVOLUME;
2568         try
2569         {
2570           OCC_CATCH_SIGNALS;
2571 #ifdef NETGEN_V5
2572           err = netgen::OCCGenerateMesh(occgeo, _ngMesh, mparams, startWith, endWith);
2573 #else
2574           err = netgen::OCCGenerateMesh(occgeo, _ngMesh, startWith, endWith, optstr);
2575 #endif
2576           if(netgen::multithread.terminate)
2577             return false;
2578
2579           if ( comment.empty() ) // do not overwrite a previos error
2580             comment << text(err);
2581         }
2582         catch (Standard_Failure& ex)
2583         {
2584           if ( comment.empty() ) // do not overwrite a previos error
2585             comment << text(ex);
2586         }
2587         catch (netgen::NgException exc)
2588         {
2589           if ( comment.empty() ) // do not overwrite a previos error
2590             comment << text(exc);
2591         }
2592       }
2593     }
2594     if (!err && mparams.secondorder > 0)
2595     {
2596       try
2597       {
2598         OCC_CATCH_SIGNALS;
2599         netgen::OCCRefinementSurfaces ref (occgeo);
2600         ref.MakeSecondOrder (*_ngMesh);
2601       }
2602       catch (Standard_Failure& ex)
2603       {
2604         if ( comment.empty() ) // do not overwrite a previos error
2605           comment << "Exception in netgen at passing to 2nd order ";
2606       }
2607       catch (netgen::NgException exc)
2608       {
2609         if ( comment.empty() ) // do not overwrite a previos error
2610           comment << exc.What();
2611       }
2612     }
2613   }
2614
2615   _ticTime = 0.98 / _progressTic;
2616
2617   int nbNod = _ngMesh->GetNP();
2618   int nbSeg = _ngMesh->GetNSeg();
2619   int nbFac = _ngMesh->GetNSE();
2620   int nbVol = _ngMesh->GetNE();
2621   bool isOK = ( !err && (_isVolume ? (nbVol > 0) : (nbFac > 0)) );
2622
2623   MESSAGE((err ? "Mesh Generation failure" : "End of Mesh Generation") <<
2624           ", nb nodes: "    << nbNod <<
2625           ", nb segments: " << nbSeg <<
2626           ", nb faces: "    << nbFac <<
2627           ", nb volumes: "  << nbVol);
2628
2629   // Feed back the SMESHDS with the generated Nodes and Elements
2630   if ( true /*isOK*/ ) // get whatever built
2631     FillSMesh( occgeo, *_ngMesh, initState, *_mesh, nodeVec, comment ); //!< 
2632
2633   SMESH_ComputeErrorPtr readErr = ReadErrors(nodeVec);
2634   if ( readErr && !readErr->myBadElements.empty() )
2635     error = readErr;
2636
2637   if ( error->IsOK() && ( !isOK || comment.size() > 0 ))
2638     error->myName = COMPERR_ALGO_FAILED;
2639   if ( !comment.empty() )
2640     error->myComment = comment;
2641
2642   // SetIsAlwaysComputed( true ) to empty sub-meshes, which
2643   // appear if the geometry contains coincident sub-shape due
2644   // to bool merge_solids = 1; in netgen/libsrc/occ/occgenmesh.cpp
2645   const int nbMaps = 2;
2646   const TopTools_IndexedMapOfShape* geoMaps[nbMaps] =
2647     { & occgeo.vmap, & occgeo.emap/*, & occgeo.fmap*/ };
2648   for ( int iMap = 0; iMap < nbMaps; ++iMap )
2649     for (int i = 1; i <= geoMaps[iMap]->Extent(); i++)
2650       if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( geoMaps[iMap]->FindKey(i)))
2651         if ( !sm->IsMeshComputed() )
2652           sm->SetIsAlwaysComputed( true );
2653
2654   // set bad compute error to subshapes of all failed sub-shapes
2655   if ( !error->IsOK() )
2656   {
2657     bool pb2D = false, pb3D = false;
2658     for (int i = 1; i <= occgeo.fmap.Extent(); i++) {
2659       int status = occgeo.facemeshstatus[i-1];
2660       if (status == 1 ) continue;
2661       if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.fmap( i ))) {
2662         SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
2663         if ( !smError || smError->IsOK() ) {
2664           if ( status == -1 )
2665             smError.reset( new SMESH_ComputeError( *error ));
2666           else
2667             smError.reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED, "Ignored" ));
2668           if ( SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
2669             smError->myName = COMPERR_WARNING;
2670         }
2671         pb2D = pb2D || smError->IsKO();
2672       }
2673     }
2674     if ( !pb2D ) // all faces are OK
2675       for (int i = 1; i <= occgeo.somap.Extent(); i++)
2676         if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( occgeo.somap( i )))
2677         {
2678           bool smComputed = nbVol && !sm->IsEmpty();
2679           if ( smComputed && internals.hasInternalVertexInSolid( sm->GetId() ))
2680           {
2681             int nbIntV = internals.getSolidsWithVertices().find( sm->GetId() )->second.size();
2682             SMESHDS_SubMesh* smDS = sm->GetSubMeshDS();
2683             smComputed = ( smDS->NbElements() > 0 || smDS->NbNodes() > nbIntV );
2684           }
2685           SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
2686           if ( !smComputed && ( !smError || smError->IsOK() ))
2687           {
2688             smError.reset( new SMESH_ComputeError( *error ));
2689             if ( nbVol && SMESH_Algo::GetMeshError( sm ) == SMESH_Algo::MEr_OK )
2690               smError->myName = COMPERR_WARNING;
2691           }
2692           pb3D = pb3D || ( smError && smError->IsKO() );
2693         }
2694     if ( !pb2D && !pb3D )
2695       err = 0; // no fatal errors, only warnings
2696   }
2697
2698   ngLib._isComputeOk = !err;
2699
2700   return !err;
2701 }
2702
2703 //=============================================================================
2704 /*!
2705  * Evaluate
2706  */
2707 //=============================================================================
2708 bool NETGENPlugin_Mesher::Evaluate(MapShapeNbElems& aResMap)
2709 {
2710   netgen::MeshingParameters& mparams = netgen::mparam;
2711
2712
2713   // -------------------------
2714   // Prepare OCC geometry
2715   // -------------------------
2716   netgen::OCCGeometry occgeo;
2717   list< SMESH_subMesh* > meshedSM[4]; // for 0-3 dimensions
2718   NETGENPlugin_Internals internals( *_mesh, _shape, _isVolume );
2719   PrepareOCCgeometry( occgeo, _shape, *_mesh, meshedSM, &internals );
2720
2721   bool tooManyElems = false;
2722   const int hugeNb = std::numeric_limits<int>::max() / 100;
2723
2724   // ----------------
2725   // evaluate 1D 
2726   // ----------------
2727   // pass 1D simple parameters to NETGEN
2728   if ( _simpleHyp )
2729   {
2730     // not to RestrictLocalH() according to curvature during MESHCONST_ANALYSE
2731     mparams.uselocalh = false;
2732     mparams.grading = 0.8; // not limitited size growth
2733
2734     if ( _simpleHyp->GetNumberOfSegments() )
2735       // nb of segments
2736       mparams.maxh = occgeo.boundingbox.Diam();
2737     else
2738       // segment length
2739       mparams.maxh = _simpleHyp->GetLocalLength();
2740   }
2741
2742   if ( mparams.maxh == 0.0 )
2743     mparams.maxh = occgeo.boundingbox.Diam();
2744   if ( _simpleHyp || ( mparams.minh == 0.0 && _fineness != NETGENPlugin_Hypothesis::UserDefined))
2745     mparams.minh = GetDefaultMinSize( _shape, mparams.maxh );
2746
2747   // let netgen create _ngMesh and calculate element size on not meshed shapes
2748   NETGENPlugin_NetgenLibWrapper ngLib;
2749   netgen::Mesh *ngMesh = NULL;
2750 #ifndef NETGEN_V5
2751   char *optstr = 0;
2752 #endif
2753   int startWith = netgen::MESHCONST_ANALYSE;
2754   int endWith   = netgen::MESHCONST_MESHEDGES;
2755 #ifdef NETGEN_V5
2756   int err = netgen::OCCGenerateMesh(occgeo, ngMesh, mparams, startWith, endWith);
2757 #else
2758   int err = netgen::OCCGenerateMesh(occgeo, ngMesh, startWith, endWith, optstr);
2759 #endif
2760
2761   if(netgen::multithread.terminate)
2762     return false;
2763
2764   ngLib.setMesh(( Ng_Mesh*) ngMesh );
2765   if (err) {
2766     if ( SMESH_subMesh* sm = _mesh->GetSubMeshContaining( _shape ))
2767       sm->GetComputeError().reset( new SMESH_ComputeError( COMPERR_ALGO_FAILED ));
2768     return false;
2769   }
2770   if ( _simpleHyp )
2771   {
2772     // Pass 1D simple parameters to NETGEN
2773     // --------------------------------
2774     int      nbSeg = _simpleHyp->GetNumberOfSegments();
2775     double segSize = _simpleHyp->GetLocalLength();
2776     for ( int iE = 1; iE <= occgeo.emap.Extent(); ++iE )
2777     {
2778       const TopoDS_Edge& e = TopoDS::Edge( occgeo.emap(iE));
2779       if ( nbSeg )
2780         segSize = SMESH_Algo::EdgeLength( e ) / ( nbSeg - 0.4 );
2781       setLocalSize( e, segSize, *ngMesh );
2782     }
2783   }
2784   else // if ( ! _simpleHyp )
2785   {
2786     // Local size on vertices and edges
2787     // --------------------------------
2788     for(std::map<int,double>::const_iterator it=EdgeId2LocalSize.begin(); it!=EdgeId2LocalSize.end(); it++)
2789     {
2790       int key = (*it).first;
2791       double hi = (*it).second;
2792       const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2793       const TopoDS_Edge& e = TopoDS::Edge(shape);
2794       setLocalSize( e, hi, *ngMesh );
2795     }
2796     for(std::map<int,double>::const_iterator it=VertexId2LocalSize.begin(); it!=VertexId2LocalSize.end(); it++)
2797     {
2798       int key = (*it).first;
2799       double hi = (*it).second;
2800       const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2801       const TopoDS_Vertex& v = TopoDS::Vertex(shape);
2802       gp_Pnt p = BRep_Tool::Pnt(v);
2803       NETGENPlugin_Mesher::RestrictLocalSize( *ngMesh, p.XYZ(), hi );
2804     }
2805     for(map<int,double>::const_iterator it=FaceId2LocalSize.begin();
2806         it!=FaceId2LocalSize.end(); it++)
2807     {
2808       int key = (*it).first;
2809       double val = (*it).second;
2810       const TopoDS_Shape& shape = ShapesWithLocalSize.FindKey(key);
2811       int faceNgID = occgeo.fmap.FindIndex(shape);
2812       occgeo.SetFaceMaxH(faceNgID, val);
2813       for ( TopExp_Explorer edgeExp( shape, TopAbs_EDGE ); edgeExp.More(); edgeExp.Next() )
2814         setLocalSize( TopoDS::Edge( edgeExp.Current() ), val, *ngMesh );
2815     }
2816   }
2817   // calculate total nb of segments and length of edges
2818   double fullLen = 0.0;
2819   int fullNbSeg = 0;
2820   int entity = mparams.secondorder > 0 ? SMDSEntity_Quad_Edge : SMDSEntity_Edge;
2821   TopTools_DataMapOfShapeInteger Edge2NbSeg;
2822   for (TopExp_Explorer exp(_shape, TopAbs_EDGE); exp.More(); exp.Next())
2823   {
2824     TopoDS_Edge E = TopoDS::Edge( exp.Current() );
2825     if( !Edge2NbSeg.Bind(E,0) )
2826       continue;
2827
2828     double aLen = SMESH_Algo::EdgeLength(E);
2829     fullLen += aLen;
2830
2831     vector<int>& aVec = aResMap[_mesh->GetSubMesh(E)];
2832     if ( aVec.empty() )
2833       aVec.resize( SMDSEntity_Last, 0);
2834     else
2835       fullNbSeg += aVec[ entity ];
2836   }
2837
2838   // store nb of segments computed by Netgen
2839   NCollection_Map<Link> linkMap;
2840   for (int i = 1; i <= ngMesh->GetNSeg(); ++i )
2841   {
2842     const netgen::Segment& seg = ngMesh->LineSegment(i);
2843     Link link(seg[0], seg[1]);
2844     if ( !linkMap.Add( link )) continue;
2845     int aGeomEdgeInd = seg.epgeominfo[0].edgenr;
2846     if (aGeomEdgeInd > 0 && aGeomEdgeInd <= occgeo.emap.Extent())
2847     {
2848       vector<int>& aVec = aResMap[_mesh->GetSubMesh(occgeo.emap(aGeomEdgeInd))];
2849       aVec[ entity ]++;
2850     }
2851   }
2852   // store nb of nodes on edges computed by Netgen
2853   TopTools_DataMapIteratorOfDataMapOfShapeInteger Edge2NbSegIt(Edge2NbSeg);
2854   for (; Edge2NbSegIt.More(); Edge2NbSegIt.Next())
2855   {
2856     vector<int>& aVec = aResMap[_mesh->GetSubMesh(Edge2NbSegIt.Key())];
2857     if ( aVec[ entity ] > 1 && aVec[ SMDSEntity_Node ] == 0 )
2858       aVec[SMDSEntity_Node] = mparams.secondorder > 0  ? 2*aVec[ entity ]-1 : aVec[ entity ]-1;
2859
2860     fullNbSeg += aVec[ entity ];
2861     Edge2NbSeg( Edge2NbSegIt.Key() ) = aVec[ entity ];
2862   }
2863   if ( fullNbSeg == 0 )
2864     return false;
2865
2866   // ----------------
2867   // evaluate 2D 
2868   // ----------------
2869   if ( _simpleHyp ) {
2870     if ( double area = _simpleHyp->GetMaxElementArea() ) {
2871       // face area
2872       mparams.maxh = sqrt(2. * area/sqrt(3.0));
2873       mparams.grading = 0.4; // moderate size growth
2874     }
2875     else {
2876       // length from edges
2877       mparams.maxh = fullLen/fullNbSeg;
2878       mparams.grading = 0.2; // slow size growth
2879     }
2880   }
2881   mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2882   mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
2883
2884   for (TopExp_Explorer exp(_shape, TopAbs_FACE); exp.More(); exp.Next())
2885   {
2886     TopoDS_Face F = TopoDS::Face( exp.Current() );
2887     SMESH_subMesh *sm = _mesh->GetSubMesh(F);
2888     GProp_GProps G;
2889     BRepGProp::SurfaceProperties(F,G);
2890     double anArea = G.Mass();
2891     tooManyElems = tooManyElems || ( anArea/hugeNb > mparams.maxh*mparams.maxh );
2892     int nb1d = 0;
2893     if ( !tooManyElems )
2894     {
2895       TopTools_MapOfShape egdes;
2896       for (TopExp_Explorer exp1(F,TopAbs_EDGE); exp1.More(); exp1.Next())
2897         if ( egdes.Add( exp1.Current() ))
2898           nb1d += Edge2NbSeg.Find(exp1.Current());
2899     }
2900     int nbFaces = tooManyElems ? hugeNb : int( 4*anArea / (mparams.maxh*mparams.maxh*sqrt(3.)));
2901     int nbNodes = tooManyElems ? hugeNb : (( nbFaces*3 - (nb1d-1)*2 ) / 6 + 1 );
2902
2903     vector<int> aVec(SMDSEntity_Last, 0);
2904     if( mparams.secondorder > 0 ) {
2905       int nb1d_in = (nbFaces*3 - nb1d) / 2;
2906       aVec[SMDSEntity_Node] = nbNodes + nb1d_in;
2907       aVec[SMDSEntity_Quad_Triangle] = nbFaces;
2908     }
2909     else {
2910       aVec[SMDSEntity_Node] = Max ( nbNodes, 0  );
2911       aVec[SMDSEntity_Triangle] = nbFaces;
2912     }
2913     aResMap[sm].swap(aVec);
2914   }
2915
2916   // ----------------
2917   // evaluate 3D
2918   // ----------------
2919   if(_isVolume) {
2920     // pass 3D simple parameters to NETGEN
2921     const NETGENPlugin_SimpleHypothesis_3D* simple3d =
2922       dynamic_cast< const NETGENPlugin_SimpleHypothesis_3D* > ( _simpleHyp );
2923     if ( simple3d ) {
2924       if ( double vol = simple3d->GetMaxElementVolume() ) {
2925         // max volume
2926         mparams.maxh = pow( 72, 1/6. ) * pow( vol, 1/3. );
2927         mparams.maxh = min( mparams.maxh, occgeo.boundingbox.Diam()/2 );
2928       }
2929       else {
2930         // using previous length from faces
2931       }
2932       mparams.grading = 0.4;
2933       mparams.maxh = min( mparams.maxh, fullLen/fullNbSeg * (1. + mparams.grading));
2934     }
2935     GProp_GProps G;
2936     BRepGProp::VolumeProperties(_shape,G);
2937     double aVolume = G.Mass();
2938     double tetrVol = 0.1179*mparams.maxh*mparams.maxh*mparams.maxh;
2939     tooManyElems = tooManyElems || ( aVolume/hugeNb > tetrVol );
2940     int nbVols = tooManyElems ? hugeNb : int(aVolume/tetrVol);
2941     int nb1d_in = int(( nbVols*6 - fullNbSeg ) / 6 );
2942     vector<int> aVec(SMDSEntity_Last, 0 );
2943     if ( tooManyElems ) // avoid FPE
2944     {
2945       aVec[SMDSEntity_Node] = hugeNb;
2946       aVec[ mparams.secondorder > 0 ? SMDSEntity_Quad_Tetra : SMDSEntity_Tetra] = hugeNb;
2947     }
2948     else
2949     {
2950       if( mparams.secondorder > 0 ) {
2951         aVec[SMDSEntity_Node] = nb1d_in/3 + 1 + nb1d_in;
2952         aVec[SMDSEntity_Quad_Tetra] = nbVols;
2953       }
2954       else {
2955         aVec[SMDSEntity_Node] = nb1d_in/3 + 1;
2956         aVec[SMDSEntity_Tetra] = nbVols;
2957       }
2958     }
2959     SMESH_subMesh *sm = _mesh->GetSubMesh(_shape);
2960     aResMap[sm].swap(aVec);
2961   }
2962
2963   return true;
2964 }
2965
2966 double NETGENPlugin_Mesher::GetProgress(const SMESH_Algo* holder,
2967                                         const int *       algoProgressTic,
2968                                         const double *    algoProgress) const
2969 {
2970   ((int&) _progressTic ) = *algoProgressTic + 1;
2971
2972   if ( !_occgeom ) return 0;
2973
2974   double progress = -1;
2975   if ( !_isVolume )
2976   {
2977     if ( _ticTime < 0 && netgen::multithread.task[0] == 'O'/*Optimizing surface*/ )
2978     {
2979       ((double&) _ticTime ) = edgeFaceMeshingTime / _totalTime / _progressTic;
2980     }
2981     else if ( !_optimize /*&& _occgeom->fmap.Extent() > 1*/ )
2982     {
2983       int doneShapeIndex = -1;
2984       while ( doneShapeIndex+1 < _occgeom->facemeshstatus.Size() &&
2985               _occgeom->facemeshstatus[ doneShapeIndex+1 ])
2986         doneShapeIndex++;
2987       if ( doneShapeIndex+1 != _curShapeIndex )
2988       {
2989         ((int&) _curShapeIndex) = doneShapeIndex+1;
2990         double    doneShapeRate = _curShapeIndex / double( _occgeom->fmap.Extent() );
2991         double         doneTime = edgeMeshingTime + doneShapeRate * faceMeshingTime;
2992         ((double&)    _ticTime) = doneTime / _totalTime / _progressTic;
2993         // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
2994         //      << " " << doneTime / _totalTime / _progressTic << endl;
2995       }
2996     }
2997   }
2998   else if ( !_optimize && _occgeom->somap.Extent() > 1 )
2999   {
3000     int curShapeIndex = _curShapeIndex;
3001     if ( _ngMesh->GetNE() > 0 )
3002     {
3003       netgen::Element el = (*_ngMesh)[netgen::ElementIndex( _ngMesh->GetNE()-1 )];
3004       curShapeIndex = el.GetIndex();
3005     }
3006     if ( curShapeIndex != _curShapeIndex )
3007     {
3008       ((int&) _curShapeIndex) = curShapeIndex;
3009       double    doneShapeRate = _curShapeIndex / double( _occgeom->somap.Extent() );
3010       double         doneTime = edgeFaceMeshingTime + doneShapeRate * voluMeshingTime;
3011       ((double&)    _ticTime) = doneTime / _totalTime / _progressTic;
3012       // cout << "shape " << _curShapeIndex << " _ticTime " << _ticTime
3013       //      << " " << doneTime / _totalTime / _progressTic << endl;
3014     }
3015   }
3016   if ( _ticTime > 0 )
3017     progress  = Max( *algoProgressTic * _ticTime, *algoProgress );
3018   if ( progress > 0 )
3019   {
3020     ((int&) *algoProgressTic )++;
3021     ((double&) *algoProgress) = progress;
3022   }
3023   //cout << progress << " "  << *algoProgressTic << " " << netgen::multithread.task << " "<< _ticTime << endl;
3024
3025   return Min( progress, 0.99 );
3026 }
3027
3028 //================================================================================
3029 /*!
3030  * \brief Remove "test.out" and "problemfaces" files in current directory
3031  */
3032 //================================================================================
3033
3034 void NETGENPlugin_Mesher::RemoveTmpFiles()
3035 {
3036   bool rm =  SMESH_File("test.out").remove() ;
3037 #ifndef WIN32
3038   if (rm && netgen::testout)
3039   {
3040     delete netgen::testout;
3041     netgen::testout = 0;
3042   }
3043 #endif
3044   SMESH_File("problemfaces").remove();
3045   SMESH_File("occmesh.rep").remove();
3046 }
3047
3048 //================================================================================
3049 /*!
3050  * \brief Read mesh entities preventing successful computation from "test.out" file
3051  */
3052 //================================================================================
3053
3054 SMESH_ComputeErrorPtr
3055 NETGENPlugin_Mesher::ReadErrors(const vector<const SMDS_MeshNode* >& nodeVec)
3056 {
3057   SMESH_ComputeErrorPtr err = SMESH_ComputeError::New
3058     (COMPERR_BAD_INPUT_MESH, "Some edges multiple times in surface mesh");
3059   SMESH_File file("test.out");
3060   vector<int> two(2);
3061   const char* badEdgeStr = " multiple times in surface mesh";
3062   const int   badEdgeStrLen = strlen( badEdgeStr );
3063   while( !file.eof() )
3064   {
3065     if ( strncmp( file, "Edge ", 5 ) == 0 &&
3066          file.getInts( two ) &&
3067          strncmp( file, badEdgeStr, badEdgeStrLen ) == 0 &&
3068          two[0] < nodeVec.size() && two[1] < nodeVec.size())
3069     {
3070       err->myBadElements.push_back( new SMDS_LinearEdge( nodeVec[ two[0]], nodeVec[ two[1]] ));
3071       file += badEdgeStrLen;
3072     }
3073     else if ( strncmp( file, "Intersecting: ", 14 ) == 0 )
3074     {
3075 // Intersecting: 
3076 // openelement 18 with open element 126
3077 // 41  36  38  
3078 // 69  70  72
3079       vector<int> three1(3), three2(3);
3080       file.getLine();
3081       const char* pos = file;
3082       bool ok = ( strncmp( file, "openelement ", 12 ) == 0 );
3083       ok = ok && file.getInts( two );
3084       ok = ok && file.getInts( three1 );
3085       ok = ok && file.getInts( three2 );
3086       for ( int i = 0; ok && i < 3; ++i )
3087         ok = ( three1[i] < nodeVec.size() && nodeVec[ three1[i]]);
3088       for ( int i = 0; ok && i < 3; ++i ) 
3089         ok = ( three2[i] < nodeVec.size() && nodeVec[ three2[i]]);
3090       if ( ok )
3091       {
3092         err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three1[0]],
3093                                                             nodeVec[ three1[1]],
3094                                                             nodeVec[ three1[2]]));
3095         err->myBadElements.push_back( new SMDS_FaceOfNodes( nodeVec[ three2[0]],
3096                                                             nodeVec[ three2[1]],
3097                                                             nodeVec[ three2[2]]));
3098         err->myComment = "Intersecting triangles";
3099       }
3100       else
3101       {
3102         file.setPos( pos );
3103       }
3104     }
3105     else
3106     {
3107       ++file;
3108     }
3109   }
3110   return err;
3111 }
3112
3113 //================================================================================
3114 /*!
3115  * \brief Write a python script creating an equivalent SALOME mesh.
3116  * This is useful to see what mesh is passed as input for the next step of mesh
3117  * generation (of mesh of higher dimension)
3118  */
3119 //================================================================================
3120
3121 void NETGENPlugin_Mesher::toPython( const netgen::Mesh* ngMesh,
3122                                     const std::string&  pyFile)
3123 {
3124   ofstream outfile(pyFile.c_str(), ios::out);
3125   if ( !outfile ) return;
3126
3127   outfile << "import smesh, SMESH" << endl
3128           << "mesh = smesh.Mesh()" << endl << endl;
3129
3130   using namespace netgen;
3131   PointIndex pi;
3132   for (pi = PointIndex::BASE; 
3133        pi < ngMesh->GetNP()+PointIndex::BASE; pi++)
3134   {
3135     outfile << "mesh.AddNode( ";
3136     outfile << (*ngMesh)[pi](0) << ", ";
3137     outfile << (*ngMesh)[pi](1) << ", ";
3138     outfile << (*ngMesh)[pi](2) << ") ## "<< pi << endl;
3139   }
3140
3141   int nbDom = ngMesh->GetNDomains();
3142   for ( int i = 0; i < nbDom; ++i )
3143     outfile<< "grp" << i+1 << " = mesh.CreateEmptyGroup( SMESH.FACE, 'domain"<< i+1 << "')"<< endl;
3144
3145   SurfaceElementIndex sei;
3146   for (sei = 0; sei < ngMesh->GetNSE(); sei++)
3147   {
3148     outfile << "mesh.AddFace([ ";
3149     Element2d sel = (*ngMesh)[sei];
3150     for (int j = 0; j < sel.GetNP(); j++)
3151       outfile << sel[j] << ( j+1 < sel.GetNP() ? ", " : " ])");
3152     if ( sel.IsDeleted() ) outfile << " ## IsDeleted ";
3153     outfile << endl;
3154
3155     if ((*ngMesh)[sei].GetIndex())
3156     {
3157       if ( int dom1 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainIn())
3158         outfile << "grp"<< dom1 <<".Add([ " << (int)sei+1 << " ])" << endl;
3159       if ( int dom2 = ngMesh->GetFaceDescriptor((*ngMesh)[sei].GetIndex ()).DomainOut())
3160         outfile << "grp"<< dom2 <<".Add([ " << (int)sei+1 << " ])" << endl;
3161     }
3162   }
3163
3164   for (ElementIndex ei = 0; ei < ngMesh->GetNE(); ei++)
3165   {
3166     Element el = (*ngMesh)[ei];
3167     outfile << "mesh.AddVolume([ ";
3168     for (int j = 0; j < el.GetNP(); j++)
3169       outfile << el[j] << ( j+1 < el.GetNP() ? ", " : " ])");
3170     outfile << endl;
3171   }
3172
3173   for (int i = 1; i <= ngMesh->GetNSeg(); i++)
3174   {
3175     const Segment & seg = ngMesh->LineSegment (i);
3176     outfile << "mesh.AddEdge([ "
3177             << seg[0] << ", "
3178             << seg[1] << " ])" << endl;
3179   }
3180   cout << "Write " << pyFile << endl;
3181 }
3182
3183 //================================================================================
3184 /*!
3185  * \brief Constructor of NETGENPlugin_ngMeshInfo
3186  */
3187 //================================================================================
3188
3189 NETGENPlugin_ngMeshInfo::NETGENPlugin_ngMeshInfo( netgen::Mesh* ngMesh):
3190   _copyOfLocalH(0)
3191 {
3192   if ( ngMesh )
3193   {
3194     _nbNodes    = ngMesh->GetNP();
3195     _nbSegments = ngMesh->GetNSeg();
3196     _nbFaces    = ngMesh->GetNSE();
3197     _nbVolumes  = ngMesh->GetNE();
3198   }
3199   else
3200   {
3201     _nbNodes = _nbSegments = _nbFaces = _nbVolumes = 0;
3202   }
3203 }
3204
3205 //================================================================================
3206 /*!
3207  * \brief Copy LocalH member from one netgen mesh to another
3208  */
3209 //================================================================================
3210
3211 void NETGENPlugin_ngMeshInfo::transferLocalH( netgen::Mesh* fromMesh,
3212                                               netgen::Mesh* toMesh )
3213 {
3214   if ( !fromMesh->LocalHFunctionGenerated() ) return;
3215   if ( !toMesh->LocalHFunctionGenerated() )
3216 #ifdef NETGEN_V5
3217     toMesh->CalcLocalH(netgen::mparam.grading);
3218 #else
3219     toMesh->CalcLocalH();
3220 #endif
3221
3222   const size_t size = sizeof( netgen::LocalH );
3223   _copyOfLocalH = new char[ size ];
3224   memcpy( (void*)_copyOfLocalH, (void*)&toMesh->LocalHFunction(), size );
3225   memcpy( (void*)&toMesh->LocalHFunction(), (void*)&fromMesh->LocalHFunction(), size );
3226 }
3227
3228 //================================================================================
3229 /*!
3230  * \brief Restore LocalH member of a netgen mesh
3231  */
3232 //================================================================================
3233
3234 void NETGENPlugin_ngMeshInfo::restoreLocalH( netgen::Mesh* toMesh )
3235 {
3236   if ( _copyOfLocalH )
3237   {
3238     const size_t size = sizeof( netgen::LocalH );
3239     memcpy( (void*)&toMesh->LocalHFunction(), (void*)_copyOfLocalH, size );
3240     delete [] _copyOfLocalH;
3241     _copyOfLocalH = 0;
3242   }
3243 }
3244
3245 //================================================================================
3246 /*!
3247  * \brief Find "internal" sub-shapes
3248  */
3249 //================================================================================
3250
3251 NETGENPlugin_Internals::NETGENPlugin_Internals( SMESH_Mesh&         mesh,
3252                                                 const TopoDS_Shape& shape,
3253                                                 bool                is3D )
3254   : _mesh( mesh ), _is3D( is3D )
3255 {
3256   SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3257
3258   TopExp_Explorer f,e;
3259   for ( f.Init( shape, TopAbs_FACE ); f.More(); f.Next() )
3260   {
3261     int faceID = meshDS->ShapeToIndex( f.Current() );
3262
3263     // find not computed internal edges
3264
3265     for ( e.Init( f.Current().Oriented(TopAbs_FORWARD), TopAbs_EDGE ); e.More(); e.Next() )
3266       if ( e.Current().Orientation() == TopAbs_INTERNAL )
3267       {
3268         SMESH_subMesh* eSM = mesh.GetSubMesh( e.Current() );
3269         if ( eSM->IsEmpty() )
3270         {
3271           _e2face.insert( make_pair( eSM->GetId(), faceID ));
3272           for ( TopoDS_Iterator v(e.Current()); v.More(); v.Next() )
3273             _e2face.insert( make_pair( meshDS->ShapeToIndex( v.Value() ), faceID ));
3274         }
3275       }
3276
3277     // find internal vertices in a face
3278     set<int> intVV; // issue 0020850 where same vertex is twice in a face
3279     for ( TopoDS_Iterator fSub( f.Current() ); fSub.More(); fSub.Next())
3280       if ( fSub.Value().ShapeType() == TopAbs_VERTEX )
3281       {
3282         int vID = meshDS->ShapeToIndex( fSub.Value() );
3283         if ( intVV.insert( vID ).second )
3284           _f2v[ faceID ].push_back( vID );
3285       }
3286
3287     if ( is3D )
3288     {
3289       // find internal faces and their subshapes where nodes are to be doubled
3290       //  to make a crack with non-sewed borders
3291
3292       if ( f.Current().Orientation() == TopAbs_INTERNAL )
3293       {
3294         _intShapes.insert( meshDS->ShapeToIndex( f.Current() ));
3295
3296         // egdes
3297         list< TopoDS_Shape > edges;
3298         for ( e.Init( f.Current(), TopAbs_EDGE ); e.More(); e.Next())
3299           if ( SMESH_MesherHelper::NbAncestors( e.Current(), mesh, TopAbs_FACE ) > 1 )
3300           {
3301             _intShapes.insert( meshDS->ShapeToIndex( e.Current() ));
3302             edges.push_back( e.Current() );
3303             // find border faces
3304             PShapeIteratorPtr fIt =
3305               SMESH_MesherHelper::GetAncestors( edges.back(),mesh,TopAbs_FACE );
3306             while ( const TopoDS_Shape* pFace = fIt->next() )
3307               if ( !pFace->IsSame( f.Current() ))
3308                 _borderFaces.insert( meshDS->ShapeToIndex( *pFace ));
3309           }
3310         // vertices
3311         // we consider vertex internal if it is shared by more than one internal edge
3312         list< TopoDS_Shape >::iterator edge = edges.begin();
3313         for ( ; edge != edges.end(); ++edge )
3314           for ( TopoDS_Iterator v( *edge ); v.More(); v.Next() )
3315           {
3316             set<int> internalEdges;
3317             PShapeIteratorPtr eIt =
3318               SMESH_MesherHelper::GetAncestors( v.Value(),mesh,TopAbs_EDGE );
3319             while ( const TopoDS_Shape* pEdge = eIt->next() )
3320             {
3321               int edgeID = meshDS->ShapeToIndex( *pEdge );
3322               if ( isInternalShape( edgeID ))
3323                 internalEdges.insert( edgeID );
3324             }
3325             if ( internalEdges.size() > 1 )
3326               _intShapes.insert( meshDS->ShapeToIndex( v.Value() ));
3327           }
3328       }
3329     }
3330   } // loop on geom faces
3331
3332   // find vertices internal in solids
3333   if ( is3D )
3334   {
3335     for ( TopExp_Explorer so(shape, TopAbs_SOLID); so.More(); so.Next())
3336     {
3337       int soID = meshDS->ShapeToIndex( so.Current() );
3338       for ( TopoDS_Iterator soSub( so.Current() ); soSub.More(); soSub.Next())
3339         if ( soSub.Value().ShapeType() == TopAbs_VERTEX )
3340           _s2v[ soID ].push_back( meshDS->ShapeToIndex( soSub.Value() ));
3341     }
3342   }
3343 }
3344
3345 //================================================================================
3346 /*!
3347  * \brief Find mesh faces on non-internal geom faces sharing internal edge
3348  * some nodes of which are to be doubled to make the second border of the "crack"
3349  */
3350 //================================================================================
3351
3352 void NETGENPlugin_Internals::findBorderElements( TIDSortedElemSet & borderElems )
3353 {
3354   if ( _intShapes.empty() ) return;
3355
3356   SMESH_Mesh& mesh = const_cast<SMESH_Mesh&>(_mesh);
3357   SMESHDS_Mesh* meshDS = mesh.GetMeshDS();
3358
3359   // loop on internal geom edges
3360   set<int>::const_iterator intShapeId = _intShapes.begin();
3361   for ( ; intShapeId != _intShapes.end(); ++intShapeId )
3362   {
3363     const TopoDS_Shape& s = meshDS->IndexToShape( *intShapeId );
3364     if ( s.ShapeType() != TopAbs_EDGE ) continue;
3365
3366     // get internal and non-internal geom faces sharing the internal edge <s>
3367     int intFace = 0;
3368     set<int>::iterator bordFace = _borderFaces.end();
3369     PShapeIteratorPtr faces = SMESH_MesherHelper::GetAncestors( s, _mesh, TopAbs_FACE );
3370     while ( const TopoDS_Shape* pFace = faces->next() )
3371     {
3372       int faceID = meshDS->ShapeToIndex( *pFace );
3373       if ( isInternalShape( faceID ))
3374         intFace = faceID;
3375       else
3376         bordFace = _borderFaces.insert( faceID ).first;
3377     }
3378     if ( bordFace == _borderFaces.end() || !intFace ) continue;
3379
3380     // get all links of mesh faces on internal geom face sharing nodes on edge <s>
3381     set< SMESH_OrientedLink > links; //!< links of faces on internal geom face
3382     list<const SMDS_MeshElement*> suspectFaces[2]; //!< mesh faces on border geom faces
3383     int nbSuspectFaces = 0;
3384     SMESHDS_SubMesh* intFaceSM = meshDS->MeshElements( intFace );
3385     if ( !intFaceSM || intFaceSM->NbElements() == 0 ) continue;
3386     SMESH_subMeshIteratorPtr smIt = mesh.GetSubMesh( s )->getDependsOnIterator(true,true);
3387     while ( smIt->more() )
3388     {
3389       SMESHDS_SubMesh* sm = smIt->next()->GetSubMeshDS();
3390       if ( !sm ) continue;
3391       SMDS_NodeIteratorPtr nIt = sm->GetNodes();
3392       while ( nIt->more() )
3393       {
3394         const SMDS_MeshNode* nOnEdge = nIt->next();
3395         SMDS_ElemIteratorPtr fIt = nOnEdge->GetInverseElementIterator(SMDSAbs_Face);
3396         while ( fIt->more() )
3397         {
3398           const SMDS_MeshElement* f = fIt->next();
3399           int nbNodes = f->NbNodes() / ( f->IsQuadratic() ? 2 : 1 );
3400           if ( intFaceSM->Contains( f ))
3401           {
3402             for ( int i = 0; i < nbNodes; ++i )
3403               links.insert( SMESH_OrientedLink( f->GetNode(i), f->GetNode((i+1)%nbNodes)));
3404           }
3405           else
3406           {
3407             int nbDblNodes = 0;
3408             for ( int i = 0; i < nbNodes; ++i )
3409               nbDblNodes += isInternalShape( f->GetNode(i)->getshapeId() );
3410             if ( nbDblNodes )
3411               suspectFaces[ nbDblNodes < 2 ].push_back( f );
3412             nbSuspectFaces++;
3413           }
3414         }
3415       }
3416     }
3417     // suspectFaces[0] having link with same orientation as mesh faces on
3418     // the internal geom face are <borderElems>. suspectFaces[1] have
3419     // only one node on edge <s>, we decide on them later (at the 2nd loop)
3420     // by links of <borderElems> found at the 1st and 2nd loops
3421     set< SMESH_OrientedLink > borderLinks;
3422     for ( int isPostponed = 0; isPostponed < 2; ++isPostponed )
3423     {
3424       list<const SMDS_MeshElement*>::iterator fIt = suspectFaces[isPostponed].begin();
3425       for ( int nbF = 0; fIt != suspectFaces[isPostponed].end(); ++fIt, ++nbF )
3426       {
3427         const SMDS_MeshElement* f = *fIt;
3428         bool isBorder = false, linkFound = false, borderLinkFound = false;
3429         list< SMESH_OrientedLink > faceLinks;
3430         int nbNodes = f->NbNodes() / ( f->IsQuadratic() ? 2 : 1 );
3431         for ( int i = 0; i < nbNodes; ++i )
3432         {
3433           SMESH_OrientedLink link( f->GetNode(i), f->GetNode((i+1)%nbNodes));
3434           faceLinks.push_back( link );
3435           if ( !linkFound )
3436           {
3437             set< SMESH_OrientedLink >::iterator foundLink = links.find( link );
3438             if ( foundLink != links.end() )
3439             {
3440               linkFound= true;
3441               isBorder = ( foundLink->_reversed == link._reversed );
3442               if ( !isBorder && !isPostponed ) break;
3443               faceLinks.pop_back();
3444             }
3445             else if ( isPostponed && !borderLinkFound )
3446             {
3447               foundLink = borderLinks.find( link );
3448               if ( foundLink != borderLinks.end() )
3449               {
3450                 borderLinkFound = true;
3451                 isBorder = ( foundLink->_reversed != link._reversed );
3452               }
3453             }
3454           }
3455         }
3456         if ( isBorder )
3457         {
3458           borderElems.insert( f );
3459           borderLinks.insert( faceLinks.begin(), faceLinks.end() );
3460         }
3461         else if ( !linkFound && !borderLinkFound )
3462         {
3463           suspectFaces[1].push_back( f );
3464           if ( nbF > 2 * nbSuspectFaces )
3465             break; // dead loop protection
3466         }
3467       }
3468     }
3469   }
3470 }
3471
3472 //================================================================================
3473 /*!
3474  * \brief put internal shapes in maps and fill in submeshes to precompute
3475  */
3476 //================================================================================
3477
3478 void NETGENPlugin_Internals::getInternalEdges( TopTools_IndexedMapOfShape& fmap,
3479                                                TopTools_IndexedMapOfShape& emap,
3480                                                TopTools_IndexedMapOfShape& vmap,
3481                                                list< SMESH_subMesh* > smToPrecompute[])
3482 {
3483   if ( !hasInternalEdges() ) return;
3484   map<int,int>::const_iterator ev_face = _e2face.begin();
3485   for ( ; ev_face != _e2face.end(); ++ev_face )
3486   {
3487     const TopoDS_Shape& ev   = _mesh.GetMeshDS()->IndexToShape( ev_face->first );
3488     const TopoDS_Shape& face = _mesh.GetMeshDS()->IndexToShape( ev_face->second );
3489
3490     ( ev.ShapeType() == TopAbs_EDGE ? emap : vmap ).Add( ev );
3491     fmap.Add( face );
3492     //cout<<"INTERNAL EDGE or VERTEX "<<ev_face->first<<" on face "<<ev_face->second<<endl;
3493
3494     smToPrecompute[ MeshDim_1D ].push_back( _mesh.GetSubMeshContaining( ev_face->first ));
3495   }
3496 }
3497
3498 //================================================================================
3499 /*!
3500  * \brief return shapes and submeshes to be meshed and already meshed boundary submeshes
3501  */
3502 //================================================================================
3503
3504 void NETGENPlugin_Internals::getInternalFaces( TopTools_IndexedMapOfShape& fmap,
3505                                                TopTools_IndexedMapOfShape& emap,
3506                                                list< SMESH_subMesh* >&     intFaceSM,
3507                                                list< SMESH_subMesh* >&     boundarySM)
3508 {
3509   if ( !hasInternalFaces() ) return;
3510
3511   // <fmap> and <emap> are for not yet meshed shapes
3512   // <intFaceSM> is for submeshes of faces
3513   // <boundarySM> is for meshed edges and vertices
3514
3515   intFaceSM.clear();
3516   boundarySM.clear();
3517
3518   set<int> shapeIDs ( _intShapes );
3519   if ( !_borderFaces.empty() )
3520     shapeIDs.insert( _borderFaces.begin(), _borderFaces.end() );
3521
3522   set<int>::const_iterator intS = shapeIDs.begin();
3523   for ( ; intS != shapeIDs.end(); ++intS )
3524   {
3525     SMESH_subMesh* sm = _mesh.GetSubMeshContaining( *intS );
3526
3527     if ( sm->GetSubShape().ShapeType() != TopAbs_FACE ) continue;
3528
3529     intFaceSM.push_back( sm );
3530
3531     // add submeshes of not computed internal faces
3532     if ( !sm->IsEmpty() ) continue;
3533
3534     SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(true,true);
3535     while ( smIt->more() )
3536     {
3537       sm = smIt->next();
3538       const TopoDS_Shape& s = sm->GetSubShape();
3539
3540       if ( sm->IsEmpty() )
3541       {
3542         // not yet meshed
3543         switch ( s.ShapeType() ) {
3544         case TopAbs_FACE: fmap.Add ( s ); break;
3545         case TopAbs_EDGE: emap.Add ( s ); break;
3546         default:;
3547         }
3548       }
3549       else
3550       {
3551         if ( s.ShapeType() != TopAbs_FACE )
3552           boundarySM.push_back( sm );
3553       }
3554     }
3555   }
3556 }
3557
3558 //================================================================================
3559 /*!
3560  * \brief Return true if given shape is to be precomputed in order to be correctly
3561  * added to netgen mesh
3562  */
3563 //================================================================================
3564
3565 bool NETGENPlugin_Internals::isShapeToPrecompute(const TopoDS_Shape& s)
3566 {
3567   int shapeID = _mesh.GetMeshDS()->ShapeToIndex( s );
3568   switch ( s.ShapeType() ) {
3569   case TopAbs_FACE  : break; //return isInternalShape( shapeID ) || isBorderFace( shapeID );
3570   case TopAbs_EDGE  : return isInternalEdge( shapeID );
3571   case TopAbs_VERTEX: break;
3572   default:;
3573   }
3574   return false;
3575 }
3576
3577 //================================================================================
3578 /*!
3579  * \brief Return SMESH
3580  */
3581 //================================================================================
3582
3583 SMESH_Mesh& NETGENPlugin_Internals::getMesh() const
3584 {
3585   return const_cast<SMESH_Mesh&>( _mesh );
3586 }
3587
3588 //================================================================================
3589 /*!
3590  * \brief Initialize netgen library
3591  */
3592 //================================================================================
3593
3594 NETGENPlugin_NetgenLibWrapper::NETGENPlugin_NetgenLibWrapper()
3595 {
3596   Ng_Init();
3597
3598   // redirect all netgen output (mycout,myerr,cout) to _outputFileName
3599   _isComputeOk    = false;
3600   _outputFileName = getOutputFileName();
3601   netgen::mycout  = new ofstream ( _outputFileName.c_str() );
3602   netgen::myerr   = netgen::mycout;
3603   _coutBuffer     = std::cout.rdbuf();
3604 #ifdef _DEBUG_
3605   cout << "NOTE: netgen output is redirected to file " << _outputFileName << endl;
3606 #endif
3607   std::cout.rdbuf( netgen::mycout->rdbuf() );
3608
3609   _ngMesh = Ng_NewMesh();
3610 }
3611
3612 //================================================================================
3613 /*!
3614  * \brief Finish using netgen library
3615  */
3616 //================================================================================
3617
3618 NETGENPlugin_NetgenLibWrapper::~NETGENPlugin_NetgenLibWrapper()
3619 {
3620   Ng_DeleteMesh( _ngMesh );
3621   Ng_Exit();
3622   NETGENPlugin_Mesher::RemoveTmpFiles();
3623   std::cout.rdbuf( _coutBuffer );
3624 #ifdef _DEBUG_
3625   if( _isComputeOk )
3626 #endif
3627     removeOutputFile();
3628 }
3629
3630 //================================================================================
3631 /*!
3632  * \brief Set netgen mesh to delete at destruction
3633  */
3634 //================================================================================
3635
3636 void NETGENPlugin_NetgenLibWrapper::setMesh( Ng_Mesh* mesh )
3637 {
3638   if ( _ngMesh )
3639     Ng_DeleteMesh( _ngMesh );
3640   _ngMesh = mesh;
3641 }
3642
3643 //================================================================================
3644 /*!
3645  * \brief Return a unique file name
3646  */
3647 //================================================================================
3648
3649 std::string NETGENPlugin_NetgenLibWrapper::getOutputFileName()
3650 {
3651   std::string aTmpDir = SALOMEDS_Tool::GetTmpDir();
3652
3653   TCollection_AsciiString aGenericName = (char*)aTmpDir.c_str();
3654   aGenericName += "NETGEN_";
3655 #ifndef WIN32
3656   aGenericName += getpid();
3657 #else
3658   aGenericName += _getpid();
3659 #endif
3660   aGenericName += "_";
3661   aGenericName += Abs((Standard_Integer)(long) aGenericName.ToCString());
3662   aGenericName += ".out";
3663
3664   return aGenericName.ToCString();
3665 }
3666
3667 //================================================================================
3668 /*!
3669  * \brief Remove file with netgen output
3670  */
3671 //================================================================================
3672
3673 void NETGENPlugin_NetgenLibWrapper::removeOutputFile()
3674 {
3675   string tmpDir = SALOMEDS_Tool::GetDirFromPath( _outputFileName );
3676   SALOMEDS::ListOfFileNames_var aFiles = new SALOMEDS::ListOfFileNames;
3677   aFiles->length(1);
3678   std::string aFileName = SALOMEDS_Tool::GetNameFromPath( _outputFileName ) + ".out";
3679   aFiles[0] = aFileName.c_str();
3680   if ( netgen::mycout)
3681   {
3682     delete netgen::mycout;
3683     netgen::mycout = 0;
3684     netgen::myerr = 0;
3685   }
3686   
3687   SALOMEDS_Tool::RemoveTemporaryFiles( tmpDir.c_str(), aFiles.in(), true );
3688 #ifdef _DEBUG_
3689   //cout << "NOTE: netgen output log was REMOVED       " << _outputFileName << endl;
3690 #endif
3691 }