]> SALOME platform Git repositories - plugins/hexoticplugin.git/blob - src/HexoticPlugin/HexoticPlugin_Hexotic.cxx
Salome HOME
588d40888578230d10d280936350189fad01c9b2
[plugins/hexoticplugin.git] / src / HexoticPlugin / HexoticPlugin_Hexotic.cxx
1 // Copyright (C) 2007-2014  CEA/DEN, EDF R&D
2 //
3 // This library is free software; you can redistribute it and/or
4 // modify it under the terms of the GNU Lesser General Public
5 // License as published by the Free Software Foundation; either
6 // version 2.1 of the License, or (at your option) any later version.
7 //
8 // This library is distributed in the hope that it will be useful,
9 // but WITHOUT ANY WARRANTY; without even the implied warranty of
10 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
11 // Lesser General Public License for more details.
12 //
13 // You should have received a copy of the GNU Lesser General Public
14 // License along with this library; if not, write to the Free Software
15 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307 USA
16 //
17 // See http://www.salome-platform.org/ or email : webmaster.salome@opencascade.com
18 //
19
20 // ---
21 // File   : HexoticPlugin_Hexotic.cxx
22 // Author : Lioka RAZAFINDRAZAKA (CEA)
23 // ---
24 //
25 #include "HexoticPlugin_Hexotic.hxx"
26 #include "HexoticPlugin_Hypothesis.hxx"
27
28 #include "utilities.h"
29
30 #ifndef WIN32
31 #include <sys/sysinfo.h>
32 #else
33 #include <errno.h>
34 #include <process.h>
35 #endif
36
37 #ifdef _DEBUG_
38 #define DUMP(txt) \
39 //  cout << txt
40 #else
41 #define DUMP(txt)
42 #endif
43
44 #include <SMESHDS_Mesh.hxx>
45 #include <SMESHDS_GroupBase.hxx>
46 #include <SMESH_ComputeError.hxx>
47 #include <SMESH_File.hxx>
48 #include <SMESH_Gen.hxx>
49 #include <SMESH_HypoFilter.hxx>
50 #include <SMESH_MesherHelper.hxx>
51 #include <SMESH_subMesh.hxx>
52
53 #include <list>
54 #include <cstdlib>
55 #include <iostream>
56
57 #include <Standard_ProgramError.hxx>
58
59 #include <BRep_Tool.hxx>
60 #include <BRepBndLib.hxx>
61 #include <BRepClass3d_SolidClassifier.hxx>
62 #include <BRepBuilderAPI_MakeVertex.hxx>
63 #include <BRepExtrema_DistShapeShape.hxx>
64 #include <OSD_File.hxx>
65 #include <Precision.hxx>
66 #include <TopExp_Explorer.hxx>
67 #include <TopTools_MapOfShape.hxx>
68 #include <TopoDS.hxx>
69 #include <gp_Pnt.hxx>
70 #include <gp_Lin.hxx>
71 #include <GCPnts_UniformAbscissa.hxx>
72 #include <IntCurvesFace_Intersector.hxx>
73 #include <BRepMesh_IncrementalMesh.hxx>
74 #include <Poly_Triangulation.hxx>
75 #include <gp_Ax3.hxx>
76 #include <GeomAdaptor_Curve.hxx>
77
78 #include <Basics_Utils.hxx>
79 #include "GEOMImpl_Types.hxx"
80 #include "GEOM_wrap.hxx"
81
82 static void removeFile( const TCollection_AsciiString& fileName )
83 {
84   try {
85     OSD_File( fileName ).Remove();
86   }
87   catch ( Standard_ProgramError ) {
88     MESSAGE("Can't remove file: " << fileName.ToCString() << " ; file does not exist or permission denied");
89   }
90 }
91
92 //=============================================================================
93 /*!
94  *  
95  */
96 //=============================================================================
97
98 HexoticPlugin_Hexotic::HexoticPlugin_Hexotic(int hypId, int studyId, SMESH_Gen* gen)
99   : SMESH_3D_Algo(hypId, studyId, gen)
100 {
101   MESSAGE("HexoticPlugin_Hexotic::HexoticPlugin_Hexotic");
102   _name = "Hexotic_3D";
103   _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);// 1 bit /shape type
104 //   _onlyUnaryInput = false;
105   _requireShape = false;
106   _iShape=0;
107   _nbShape=0;
108   _hexoticFilesKept=false;
109   _compatibleHypothesis.push_back("Hexotic_Parameters");
110 #ifdef WITH_BLSURFPLUGIN
111   _blsurfHypo = NULL;
112 #endif
113   _compute_canceled = false;
114   
115   // Copy of what is done in BLSURFPLugin TODO : share the code
116   smeshGen_i = SMESH_Gen_i::GetSMESHGen();
117   CORBA::Object_var anObject = smeshGen_i->GetNS()->Resolve("/myStudyManager");
118   SALOMEDS::StudyManager_var aStudyMgr = SALOMEDS::StudyManager::_narrow(anObject);
119   
120   myStudy = NULL;
121   myStudy = aStudyMgr->GetStudyByID(_studyId);
122   if (myStudy)
123     MESSAGE("myStudy->StudyId() = " << myStudy->StudyId());
124 }
125
126 //=============================================================================
127 /*!
128  *  
129  */
130 //=============================================================================
131
132 HexoticPlugin_Hexotic::~HexoticPlugin_Hexotic()
133 {
134   MESSAGE("HexoticPlugin_Hexotic::~HexoticPlugin_Hexotic");
135 }
136
137
138 #ifdef WITH_BLSURFPLUGIN
139 bool HexoticPlugin_Hexotic::CheckBLSURFHypothesis( SMESH_Mesh&         aMesh,
140                                                    const TopoDS_Shape& aShape )
141 {
142   // MESSAGE("HexoticPlugin_Hexotic::CheckBLSURFHypothesis");
143   _blsurfHypo = NULL;
144
145   std::list<const SMESHDS_Hypothesis*>::const_iterator itl;
146   const SMESHDS_Hypothesis* theHyp;
147
148   // If a BLSURF hypothesis is applied, get it
149   SMESH_HypoFilter blsurfFilter;
150   blsurfFilter.Init( blsurfFilter.HasName( "BLSURF_Parameters" ));
151   std::list<const SMESHDS_Hypothesis *> appliedHyps;
152   aMesh.GetHypotheses( aShape, blsurfFilter, appliedHyps, false );
153
154   if ( appliedHyps.size() > 0 ) {
155     itl = appliedHyps.begin();
156     theHyp = (*itl); // use only the first hypothesis
157     std::string hypName = theHyp->GetName();
158     if (hypName == "BLSURF_Parameters") {
159       _blsurfHypo = static_cast<const BLSURFPlugin_Hypothesis*> (theHyp);
160       ASSERT(_blsurfHypo);
161       return true;
162     }
163   }
164   return false;
165 }
166 #endif
167
168 //=============================================================================
169 /*!
170  *  
171  */
172 //=============================================================================
173
174 bool HexoticPlugin_Hexotic::CheckHypothesis( SMESH_Mesh&                          aMesh,
175                                              const TopoDS_Shape&                  aShape,
176                                              SMESH_Hypothesis::Hypothesis_Status& aStatus )
177 {
178   // MESSAGE("HexoticPlugin_Hexotic::CheckHypothesis");
179   _hypothesis = NULL;
180
181   std::list<const SMESHDS_Hypothesis*>::const_iterator itl;
182   const SMESHDS_Hypothesis* theHyp;
183
184   const std::list<const SMESHDS_Hypothesis*>& hyps = GetUsedHypothesis(aMesh, aShape, false);
185   int nbHyp = hyps.size();
186   if (!nbHyp) {
187     aStatus = SMESH_Hypothesis::HYP_OK;
188     return true;  // can work with no hypothesis
189   }
190
191   itl = hyps.begin();
192   theHyp = (*itl); // use only the first hypothesis
193
194   std::string hypName = theHyp->GetName();
195   if (hypName == "Hexotic_Parameters") {
196     _hypothesis = static_cast<const HexoticPlugin_Hypothesis*> (theHyp);
197     ASSERT(_hypothesis);
198     aStatus = SMESH_Hypothesis::HYP_OK;
199   }
200   else
201     aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
202   
203 #ifdef WITH_BLSURFPLUGIN
204   CheckBLSURFHypothesis(aMesh, aShape);
205 #endif
206   
207   return aStatus == SMESH_Hypothesis::HYP_OK;
208 }
209
210 //=======================================================================
211 //function : findShape
212 //purpose  :
213 //=======================================================================
214
215 static TopoDS_Shape findShape(SMDS_MeshNode**     t_Node,
216                               TopoDS_Shape        aShape,
217                               const TopoDS_Shape* t_Shape,
218                               double**            t_Box,
219                               const int           nShape)
220 {
221   double pntCoor[3];
222   int iShape, nbNode = 8;
223
224   for ( int i=0; i<3; i++ ) {
225     pntCoor[i] = 0;
226     for ( int j=0; j<nbNode; j++ ) {
227       if ( i == 0) pntCoor[i] += t_Node[j]->X();
228       if ( i == 1) pntCoor[i] += t_Node[j]->Y();
229       if ( i == 2) pntCoor[i] += t_Node[j]->Z();
230     }
231     pntCoor[i] /= nbNode;
232   }
233   gp_Pnt aPnt(pntCoor[0], pntCoor[1], pntCoor[2]);
234
235   if ( aShape.IsNull() ) aShape = t_Shape[0];
236   BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
237   if ( !(SC.State() == TopAbs_IN) ) {
238     aShape.Nullify();
239     for (iShape = 0; iShape < nShape && aShape.IsNull(); iShape++) {
240       if ( !( pntCoor[0] < t_Box[iShape][0] || t_Box[iShape][1] < pntCoor[0] ||
241               pntCoor[1] < t_Box[iShape][2] || t_Box[iShape][3] < pntCoor[1] ||
242               pntCoor[2] < t_Box[iShape][4] || t_Box[iShape][5] < pntCoor[2]) ) {
243         BRepClass3d_SolidClassifier SC (t_Shape[iShape], aPnt, Precision::Confusion());
244         if (SC.State() == TopAbs_IN)
245           aShape = t_Shape[iShape];
246       }
247     }
248   }
249   return aShape;
250 }
251
252 //=======================================================================
253 //function : findEdge
254 //purpose  :
255 //=======================================================================
256
257 static int findEdge(const SMDS_MeshNode* aNode,
258                     const SMESHDS_Mesh*  aMesh,
259                     const int            nEdge,
260                     const TopoDS_Shape*  t_Edge) {
261
262   TopoDS_Shape aPntShape, foundEdge;
263   TopoDS_Vertex aVertex;
264   gp_Pnt aPnt( aNode->X(), aNode->Y(), aNode->Z() );
265
266   int foundInd, ind;
267   double nearest = RealLast(), *t_Dist;
268   double epsilon = Precision::Confusion();
269
270   t_Dist = new double[ nEdge ];
271   aPntShape = BRepBuilderAPI_MakeVertex( aPnt ).Shape();
272   aVertex   = TopoDS::Vertex( aPntShape );
273
274   for ( ind=0; ind < nEdge; ind++ ) {
275     BRepExtrema_DistShapeShape aDistance ( aVertex, t_Edge[ind] );
276     t_Dist[ind] = aDistance.Value();
277     if ( t_Dist[ind] < nearest ) {
278       nearest   = t_Dist[ind];
279       foundEdge = t_Edge[ind];
280       foundInd  = ind;
281       if ( nearest < epsilon )
282         ind = nEdge;
283     }
284   }
285
286   delete [] t_Dist;
287   return aMesh->ShapeToIndex( foundEdge );
288 }
289
290 //=======================================================================
291 //function : getNbShape
292 //purpose  :
293 //=======================================================================
294
295 static int getNbShape(std::string aFile, std::string aString, int defaultValue=0) {
296   int number = defaultValue;
297   std::string aLine;
298   std::ifstream file(aFile.c_str());
299   while ( !file.eof() ) {
300     getline( file, aLine);
301     if ( aLine == aString ) {
302       getline( file, aLine);
303       std::istringstream stringFlux( aLine );
304       stringFlux >> number;
305       number = ( number + defaultValue + std::abs(number - defaultValue) ) / 2;
306       break;
307     }
308   }
309   file.close();
310   return number;
311 }
312
313 //=======================================================================
314 //function : countShape
315 //purpose  :
316 //=======================================================================
317
318 template < class Mesh, class Shape >
319 static int countShape( Mesh* mesh, Shape shape ) {
320   TopExp_Explorer expShape ( mesh->ShapeToMesh(), shape );
321   TopTools_MapOfShape mapShape;
322   int nbShape = 0;
323   for ( ; expShape.More(); expShape.Next() ) {
324     if (mapShape.Add(expShape.Current())) {
325       nbShape++;
326     }
327   }
328   return nbShape;
329 }
330
331 //=======================================================================
332 //function : getShape
333 //purpose  :
334 //=======================================================================
335
336 template < class Mesh, class Shape, class Tab >
337 void getShape(Mesh* mesh, Shape shape, Tab *t_Shape) {
338   TopExp_Explorer expShape ( mesh->ShapeToMesh(), shape );
339   TopTools_MapOfShape mapShape;
340   for ( int i=0; expShape.More(); expShape.Next() ) {
341     if (mapShape.Add(expShape.Current())) {
342       t_Shape[i] = expShape.Current();
343       i++;
344     }
345   }
346   return;
347 }
348
349 //=======================================================================
350 //function : printWarning
351 //purpose  :
352 //=======================================================================
353
354 static void printWarning(const int nbExpected, std::string aString, const int nbFound) {
355   cout << std::endl;
356   cout << "WARNING : " << nbExpected << " " << aString << " expected, Hexotic has found " << nbFound << std::endl;
357   cout << "=======" << std::endl;
358   cout << std::endl;
359   return;
360 }
361
362 //=======================================================================
363 //function : removeHexoticFiles
364 //purpose  :
365 //=======================================================================
366
367 static void removeHexoticFiles(TCollection_AsciiString file_In, TCollection_AsciiString file_Out) {
368   removeFile( file_In );
369   removeFile( file_Out );
370 }
371
372 //=======================================================================
373 //function : splitQuads
374 //purpose  : splits all quadrangles into triangles
375 //=======================================================================
376
377 static void splitQuads(SMESH_Mesh& aMesh)
378 {
379   SMESH_MeshEditor spliter( &aMesh );
380
381   TIDSortedElemSet elems;
382   SMDS_ElemIteratorPtr eIt = aMesh.GetMeshDS()->elementsIterator();
383   while( eIt->more() )
384     elems.insert( elems.end(), eIt->next() );
385   
386   spliter.QuadToTri ( elems, /*the13Diag=*/true);
387 }
388
389 //=======================================================================
390 //function : readResult
391 //purpose  : Read GMF file in case of a mesh with geometry
392 //=======================================================================
393
394 static bool readResult(std::string         theFile,
395                        HexoticPlugin_Hexotic*  theAlgo,
396                        SMESHDS_Mesh*       theMesh,
397                        const int           nbShape,
398                        const TopoDS_Shape* tabShape,
399                        double**            tabBox)
400 {
401   // ---------------------------------
402   // Optimisation of the plugin ...
403   // Retrieve the correspondance edge --> shape
404   // (which is very costly) only when user
405   // has defined at least one group of edges
406   // which should be rare for a 3d mesh !
407   // ---------------------------------
408   
409   bool retrieve_edges = false;
410   const std::set<SMESHDS_GroupBase*>& aGroups = theMesh->GetGroups();
411   set<SMESHDS_GroupBase*>::const_iterator GrIt = aGroups.begin();
412   for (; GrIt != aGroups.end(); GrIt++)
413     {
414       SMESHDS_GroupBase* aGrp = *GrIt;
415       if ( !aGrp )
416         continue;
417       if ( aGrp->GetType() == SMDSAbs_Edge )
418         {
419           retrieve_edges = true;
420           break;
421         }
422     }
423   
424   // ---------------------------------
425   // Read generated elements and nodes
426   // ---------------------------------
427
428   TopoDS_Shape aShape;
429   TopoDS_Vertex aVertex;
430   std::string token;
431   int EndOfFile = 0, nbElem = 0, nField = 9, nbRef = 0;
432   int aHexoticNodeID = 0, shapeID, hexoticShapeID;
433   const int IdShapeRef = 2;
434   int *tabID, *tabRef, *nodeAssigne;
435   bool *tabDummy, hasDummy = false;
436   double epsilon = Precision::Confusion();
437   std::map <std::string,int> mapField;
438   SMDS_MeshNode** HexoticNode;
439   TopoDS_Shape *tabCorner, *tabEdge;
440
441   const int nbDomains = countShape( theMesh, TopAbs_SHELL );
442   const int holeID = -1;
443
444   // tabID    = new int[nbShape];
445   tabID    = new int[nbDomains];
446   tabRef   = new int[nField];
447   tabDummy = new bool[nField];
448
449   for (int i=0; i<nbDomains; i++)
450     tabID[i] = 0;
451   if ( nbDomains == 1 )
452     tabID[0] = theMesh->ShapeToIndex( tabShape[0] );
453
454   mapField["MeshVersionFormatted"] = 0; tabRef[0] = 0; tabDummy[0] = false;
455   mapField["Dimension"]            = 1; tabRef[1] = 0; tabDummy[1] = false;
456   mapField["Vertices"]             = 2; tabRef[2] = 3; tabDummy[2] = true;
457   mapField["Corners"]              = 3; tabRef[3] = 1; tabDummy[3] = false;
458   mapField["Edges"]                = 4; tabRef[4] = 2; tabDummy[4] = true;
459   mapField["Ridges"]               = 5; tabRef[5] = 1; tabDummy[5] = false;
460   mapField["Quadrilaterals"]       = 6; tabRef[6] = 4; tabDummy[6] = true;
461   mapField["Hexahedra"]            = 7; tabRef[7] = 8; tabDummy[7] = true;
462   mapField["End"]                  = 8; tabRef[8] = 0; tabDummy[0] = false;
463
464   SMDS_NodeIteratorPtr itOnHexoticInputNode = theMesh->nodesIterator();
465   while ( itOnHexoticInputNode->more() )
466     theMesh->RemoveNode( itOnHexoticInputNode->next() );
467
468   int nbVertices   = getNbShape(theFile, "Vertices");
469   int nbCorners    = getNbShape(theFile, "Corners", countShape( theMesh, TopAbs_VERTEX ));
470   int nbShapeEdge  = countShape( theMesh, TopAbs_EDGE );
471
472   tabCorner   = new TopoDS_Shape[ nbCorners ];
473   tabEdge     = new TopoDS_Shape[ nbShapeEdge ];
474   nodeAssigne = new int[ nbVertices + 1 ];
475   HexoticNode = new SMDS_MeshNode*[ nbVertices + 1 ];
476
477   getShape(theMesh, TopAbs_VERTEX, tabCorner);
478   getShape(theMesh, TopAbs_EDGE,   tabEdge);
479
480   MESSAGE("Read " << theFile << " file");
481   std::ifstream fileRes(theFile.c_str());
482   ASSERT(fileRes);
483
484   while ( EndOfFile == 0  ) {
485     int dummy;
486     fileRes >> token;
487
488     if (mapField.count(token)) {
489       nField   = mapField[token];
490       nbRef    = tabRef[nField];
491       hasDummy = tabDummy[nField];
492     }
493     else {
494       nField = -1;
495       nbRef = 0;
496     }
497
498     nbElem = 0;
499     if ( nField < (mapField.size() - 1) && nField >= 0 )
500       fileRes >> nbElem;
501
502     switch (nField) {
503       case 0: { // "MeshVersionFormatted"
504         MESSAGE(token << " " << nbElem);
505         break;
506       }
507       case 1: { // "Dimension"
508         MESSAGE("Mesh dimension " << nbElem << "D");
509         break;
510       }
511       case 2: { // "Vertices"
512         MESSAGE("Read " << nbElem << " " << token);
513         int aHexoticID;
514         double *coord;
515         SMDS_MeshNode * aHexoticNode;
516
517         coord = new double[nbRef];
518         for ( int iElem = 0; iElem < nbElem; iElem++ ) {
519           if(theAlgo->computeCanceled())
520             {
521               return false;
522             }
523           aHexoticID = iElem + 1;
524           for ( int iCoord = 0; iCoord < 3; iCoord++ )
525             fileRes >> coord[ iCoord ];
526           fileRes >> dummy;
527           aHexoticNode = theMesh->AddNode(coord[0], coord[1], coord[2]);
528           HexoticNode[ aHexoticID ] = aHexoticNode;
529           nodeAssigne[ aHexoticID ] = 0;
530         }
531         delete [] coord;
532         break;
533       }
534       case 3: // "Corners"
535       case 4: // "Edges"
536       case 5: // "Ridges"
537       case 6: // "Quadrilaterals"
538       case 7: { // "Hexahedra"
539         MESSAGE("Read " << nbElem << " " << token);
540         SMDS_MeshNode** node;
541         int nodeDim, *nodeID;
542         SMDS_MeshElement * aHexoticElement = 0;
543
544         node   = new SMDS_MeshNode*[ nbRef ];
545         nodeID = new int[ nbRef ];
546         for ( int iElem = 0; iElem < nbElem; iElem++ ) {
547           if(theAlgo->computeCanceled())
548             {
549               return false;
550             }
551           for ( int iRef = 0; iRef < nbRef; iRef++ ) {
552             fileRes >> aHexoticNodeID;                          // read nbRef aHexoticNodeID
553             node[ iRef ]   = HexoticNode[ aHexoticNodeID ];
554             nodeID[ iRef ] = aHexoticNodeID;
555           }
556           if ( hasDummy )
557             fileRes >> dummy;
558           switch (nField) {
559             case 3: { // "Corners"
560               nodeDim = 1;
561               gp_Pnt HexoticPnt ( node[0]->X(), node[0]->Y(), node[0]->Z() );
562               for ( int i=0; i<nbElem; i++ ) {
563                 aVertex = TopoDS::Vertex( tabCorner[i] );
564                 gp_Pnt aPnt = BRep_Tool::Pnt( aVertex );
565                 if ( aPnt.Distance( HexoticPnt ) < epsilon )
566                   break;
567               }
568               break;
569             }
570             case 4: { // "Edges"
571               nodeDim = 2;
572               aHexoticElement = theMesh->AddEdge( node[0], node[1] );
573               int iNode = 1;
574               if ( nodeAssigne[ nodeID[0] ] == 0 || nodeAssigne[ nodeID[0] ] == 2 )
575                 iNode = 0;
576               if(retrieve_edges)
577                 shapeID = findEdge( node[iNode], theMesh, nbShapeEdge, tabEdge );
578               else
579                 shapeID = 0;
580               break;
581             }
582             case 5: { // "Ridges"
583               break;
584             }
585             case 6: { // "Quadrilaterals"
586               nodeDim = 3;
587               aHexoticElement = theMesh->AddFace( node[0], node[1], node[2], node[3] );
588               shapeID = dummy;
589               break;
590             }
591             case 7: { // "Hexahedra"
592               nodeDim = 4;
593               if ( nbDomains > 1 ) {
594                 hexoticShapeID = dummy - IdShapeRef;
595                 if ( tabID[ hexoticShapeID ] == 0 ) {
596                   aShape = findShape(node, aShape, tabShape, tabBox, nbShape);
597                   shapeID = aShape.IsNull() ? holeID : theMesh->ShapeToIndex( aShape );
598                   tabID[ hexoticShapeID ] = shapeID;
599                 }
600                 else
601                   shapeID = tabID[ hexoticShapeID ];
602                 if ( iElem == (nbElem - 1) ) {
603                   int shapeAssociated = 0;
604                   for ( int i=0; i<nbDomains; i++ ) {
605                     if (tabID[i] > 0 )
606                       shapeAssociated += 1;
607                   }
608                   if ( shapeAssociated != nbShape )
609                     printWarning(nbShape, "domains", shapeAssociated);
610                 }
611               }
612               else {
613                 shapeID = tabID[0];
614               }
615               if ( shapeID != holeID )
616                 aHexoticElement = theMesh->AddVolume( node[0], node[3], node[2], node[1], node[4], node[7], node[6], node[5] );
617               break;
618             }
619           } // switch (nField)
620
621           if ( token != "Ridges" && ( shapeID > 0 || token == "Corners")) {
622             for ( int i=0; i<nbRef; i++ ) {
623               if ( nodeAssigne[ nodeID[i] ] == 0 ) {
624                 if      ( token == "Corners" )        theMesh->SetNodeOnVertex( node[0], aVertex );
625                 else if ( token == "Edges" )          theMesh->SetNodeOnEdge( node[i], shapeID );
626                 else if ( token == "Quadrilaterals" ) theMesh->SetNodeOnFace( node[i], shapeID );
627                 else if ( token == "Hexahedra" )      theMesh->SetNodeInVolume( node[i], shapeID );
628                 nodeAssigne[ nodeID[i] ] = nodeDim;
629               }
630             }
631             if ( token != "Corners" && aHexoticElement )
632               theMesh->SetMeshElementOnShape( aHexoticElement, shapeID );
633           }
634         }
635         delete [] node;
636         delete [] nodeID;
637         break;
638       }
639       case 8: { // "End"
640         EndOfFile = 1;
641         MESSAGE("End of " << theFile << " file");
642         break;
643       }
644       default: {
645         MESSAGE("Unknown Token: " << token);
646       }
647     }
648   }
649   cout << std::endl;
650
651   // remove nodes in holes
652   if ( nbDomains > 1 )
653   {
654     SMESHDS_SubMesh* subMesh;
655     for ( int i = 1; i <= nbVertices; ++i )
656       if ( HexoticNode[i]->NbInverseElements() == 0 )
657       {
658         subMesh =  HexoticNode[i]->getshapeId() > 0 ? theMesh->MeshElements(HexoticNode[i]->getshapeId() ) : 0;
659         theMesh->RemoveFreeNode( HexoticNode[i], subMesh, /*fromGroups=*/false );
660       }
661   }
662   delete [] tabID;
663   delete [] tabRef;
664   delete [] tabDummy;
665   delete [] tabCorner;
666   delete [] tabEdge;
667   delete [] nodeAssigne;
668   delete [] HexoticNode;
669   return true;
670 }
671
672
673 //=======================================================================
674 //function : readResult
675 //purpose  : Read GMF file in case of a mesh w/o geometry
676 //=======================================================================
677
678 static bool readResult(std::string theFile,
679                        HexoticPlugin_Hexotic*  theAlgo,
680                        SMESH_MesherHelper* theHelper)
681 {
682   SMESHDS_Mesh* theMesh = theHelper->GetMeshDS();
683
684   // ---------------------------------
685   // Read generated elements and nodes
686   // ---------------------------------
687
688   std::string token;
689   const int nbField = 9;
690   int nField, EndOfFile = 0, nbElem = 0, nbRef = 0;
691   int aHexoticNodeID = 0, shapeID;
692   int tabRef[nbField], *nodeAssigne;
693   bool tabDummy[nbField], hasDummy = false;
694   std::map <std::string,int> mapField;
695   SMDS_MeshNode** HexoticNode;
696
697   mapField["MeshVersionFormatted"] = 0; tabRef[0] = 0; tabDummy[0] = false;
698   mapField["Dimension"]            = 1; tabRef[1] = 0; tabDummy[1] = false;
699   mapField["Vertices"]             = 2; tabRef[2] = 3; tabDummy[2] = true;
700   mapField["Corners"]              = 3; tabRef[3] = 1; tabDummy[3] = false;
701   mapField["Edges"]                = 4; tabRef[4] = 2; tabDummy[4] = true;
702   mapField["Ridges"]               = 5; tabRef[5] = 1; tabDummy[5] = false;
703   mapField["Quadrilaterals"]       = 6; tabRef[6] = 4; tabDummy[6] = true;
704   mapField["Hexahedra"]            = 7; tabRef[7] = 8; tabDummy[7] = true;
705   mapField["End"]                  = 8; tabRef[8] = 0; tabDummy[8] = false;
706
707   {
708     // theMesh->Clear(); -- this does not remove imported mesh
709     SMDS_ElemIteratorPtr eIt = theMesh->elementsIterator();
710     while( eIt->more() )
711       theMesh->RemoveFreeElement( eIt->next(), /*sm=*/0 );
712     SMDS_NodeIteratorPtr nIt = theMesh->nodesIterator();
713     while ( nIt->more() )
714       theMesh->RemoveFreeNode( nIt->next(), /*sm=*/0 );
715   }
716
717   int nbVertices = getNbShape(theFile, "Vertices");
718   HexoticNode = new SMDS_MeshNode*[ nbVertices + 1 ];
719   nodeAssigne = new int[ nbVertices + 1 ];
720
721   MESSAGE("Read " << theFile << " file");
722   std::ifstream fileRes(theFile.c_str());
723   ASSERT(fileRes);
724
725   while ( !EndOfFile  )
726   {
727     int dummy;
728     fileRes >> token;
729
730     if (mapField.count(token)) {
731       nField   = mapField[token];
732       nbRef    = tabRef[nField];
733       hasDummy = tabDummy[nField];
734     }
735     else {
736       nField = -1;
737       nbRef = 0;
738     }
739
740     nbElem = 0;
741     if ( nField < (mapField.size() - 1) && nField >= 0 )
742       fileRes >> nbElem;
743
744     switch (nField) {
745     case 0: { // "MeshVersionFormatted"
746       MESSAGE(token << " " << nbElem);
747       break;
748     }
749     case 1: { // "Dimension"
750       MESSAGE("Mesh dimension " << nbElem << "D");
751       break;
752     }
753     case 2: { // "Vertices"
754       MESSAGE("Read " << nbElem << " " << token);
755       int aHexoticID;
756       double coord[3];
757       SMDS_MeshNode * aHexoticNode;
758
759       for ( int iElem = 0; iElem < nbElem; iElem++ ) {
760         if(theAlgo->computeCanceled())
761           {
762             return false;
763           }
764         aHexoticID = iElem + 1;
765         for ( int iCoord = 0; iCoord < 3; iCoord++ )
766           fileRes >> coord[ iCoord ];
767         fileRes >> dummy;
768         aHexoticNode = theMesh->AddNode(coord[0], coord[1], coord[2]);
769         HexoticNode[ aHexoticID ] = aHexoticNode;
770         nodeAssigne[ aHexoticID ] = 0;
771       }
772       break;
773     }
774     case 3: // "Corners"
775     case 4: // "Edges"
776     case 5: // "Ridges"
777     case 6: // "Quadrilaterals"
778     case 7: { // "Hexahedra"
779       MESSAGE("Read " << nbElem << " " << token);
780       std::vector< SMDS_MeshNode* > node( nbRef );
781       std::vector< int >          nodeID( nbRef );
782
783       for ( int iElem = 0; iElem < nbElem; iElem++ )
784       {
785         if(theAlgo->computeCanceled())
786           {
787             return false;
788           }
789         for ( int iRef = 0; iRef < nbRef; iRef++ )
790         {
791           fileRes >> aHexoticNodeID;                          // read nbRef aHexoticNodeID
792           node  [ iRef ] = HexoticNode[ aHexoticNodeID ];
793           nodeID[ iRef ] = aHexoticNodeID;
794         }
795         if ( hasDummy )
796           fileRes >> dummy;
797         switch (nField)
798         {
799         case 4: // "Edges"
800           theHelper->AddEdge( node[0], node[1] ); break;
801         case 6:  // "Quadrilaterals"
802           theMesh->AddFace( node[0], node[1], node[2], node[3] ); break;
803         case 7: // "Hexahedra"
804           theHelper->AddVolume( node[0], node[3], node[2], node[1],
805                                 node[4], node[7], node[6], node[5] ); break;
806         default: continue;
807         }
808         if ( nField == 6 )
809           for ( int iRef = 0; iRef < nbRef; iRef++ )
810             nodeAssigne[ nodeID[ iRef ]] = 1;
811       }
812       break;
813     }
814     case 8: { // "End"
815       EndOfFile = 1;
816       MESSAGE("End of " << theFile << " file");
817       break;
818     }
819     default: {
820       MESSAGE("Unknown Token: " << token);
821     }
822     }
823   }
824   cout << std::endl;
825
826   shapeID = theHelper->GetSubShapeID();
827   for ( int i = 0; i < nbVertices; ++i )
828     if ( !nodeAssigne[ i+1 ])
829       theMesh->SetNodeInVolume( HexoticNode[ i+1 ], shapeID );
830
831   delete [] HexoticNode;
832   delete [] nodeAssigne;
833   return true;
834 }
835
836 //=============================================================================
837 /*!
838  * Pass parameters to Hexotic
839  */
840 //=============================================================================
841
842 void HexoticPlugin_Hexotic::SetParameters(const HexoticPlugin_Hypothesis* hyp) {
843
844   MESSAGE("HexoticPlugin_Hexotic::SetParameters");
845   if (hyp) {
846     _hexesMinLevel = hyp->GetHexesMinLevel();
847     _hexesMaxLevel = hyp->GetHexesMaxLevel();
848     _hexesMinSize = hyp->GetMinSize();
849     _hexesMaxSize = hyp->GetMaxSize();
850     _hexoticIgnoreRidges = hyp->GetHexoticIgnoreRidges();
851     _hexoticInvalidElements = hyp->GetHexoticInvalidElements();
852     _hexoticSharpAngleThreshold = hyp->GetHexoticSharpAngleThreshold();
853     _hexoticNbProc = hyp->GetHexoticNbProc();
854     _hexoticWorkingDirectory = hyp->GetHexoticWorkingDirectory();
855     _hexoticVerbosity = hyp->GetHexoticVerbosity();
856     _hexoticMaxMemory = hyp->GetHexoticMaxMemory();
857     _hexoticSdMode = hyp->GetHexoticSdMode();
858     _sizeMaps = hyp->GetSizeMaps();
859   }
860   else {
861     cout << std::endl;
862     cout << "WARNING : The Hexotic default parameters are taken into account" << std::endl;
863     cout << "=======" << std::endl;
864     _hexesMinLevel = hyp->GetDefaultHexesMinLevel();
865     _hexesMaxLevel = hyp->GetDefaultHexesMaxLevel();
866     _hexesMinSize = hyp->GetDefaultMinSize();
867     _hexesMaxSize = hyp->GetDefaultMaxSize();
868     _hexoticIgnoreRidges = hyp->GetDefaultHexoticIgnoreRidges();
869     _hexoticInvalidElements = hyp->GetDefaultHexoticInvalidElements();
870     _hexoticSharpAngleThreshold = hyp->GetDefaultHexoticSharpAngleThreshold();
871     _hexoticNbProc = hyp->GetDefaultHexoticNbProc();
872     _hexoticWorkingDirectory = hyp->GetDefaultHexoticWorkingDirectory();
873     _hexoticVerbosity = hyp->GetDefaultHexoticVerbosity();
874     _hexoticMaxMemory = hyp->GetDefaultHexoticMaxMemory();
875     _hexoticSdMode = hyp->GetDefaultHexoticSdMode();
876     _sizeMaps = hyp->GetDefaultHexoticSizeMaps();
877   }
878 }
879
880 //=======================================================================
881 //function : getTmpDir
882 //purpose  :
883 //=======================================================================
884
885 static TCollection_AsciiString getTmpDir()
886 {
887   TCollection_AsciiString aTmpDir;
888
889   char *Tmp_dir = getenv("SALOME_TMP_DIR");
890 #ifdef WIN32
891   if(Tmp_dir == NULL) {
892     Tmp_dir = getenv("TEMP");
893     if( Tmp_dir== NULL )
894       Tmp_dir = getenv("TMP");
895   }
896 #endif
897
898   if(Tmp_dir != NULL) {
899     aTmpDir = Tmp_dir;
900 #ifdef WIN32
901     if(aTmpDir.Value(aTmpDir.Length()) != '\\') aTmpDir+='\\';
902 #else
903     if(aTmpDir.Value(aTmpDir.Length()) != '/') aTmpDir+='/';
904 #endif
905   }
906   else {
907 #ifdef WIN32
908     aTmpDir = TCollection_AsciiString("C:\\");
909 #else
910     aTmpDir = TCollection_AsciiString("/tmp/");
911 #endif
912   }
913   return aTmpDir;
914 }
915
916 //=======================================================================
917 //function : getSuffix
918 //purpose  : Returns a suffix that will be unique for the current process
919 //=======================================================================
920
921 static TCollection_AsciiString getSuffix()
922 {
923   TCollection_AsciiString aSuffix = "";
924   aSuffix += "_";
925 #ifndef WIN32
926   aSuffix += getenv("USER");
927 #else
928   aSuffix += getenv("USERNAME");
929 #endif
930   aSuffix += "_";
931   aSuffix += Kernel_Utils::GetHostname().c_str();
932   aSuffix += "_";
933 #ifndef WIN32
934   aSuffix += getpid();
935 #else
936   aSuffix += _getpid();
937 #endif
938
939   return aSuffix;
940 }
941
942 //================================================================================
943 /*!
944  * \brief Returns a command to run Hexotic mesher
945  */
946 //================================================================================
947
948 std::string HexoticPlugin_Hexotic::getHexoticCommand(const TCollection_AsciiString& Hexotic_In,
949                                                      const TCollection_AsciiString& Hexotic_Out,
950                                                      const TCollection_AsciiString& Hexotic_SizeMap_Prefix) const
951 {
952   cout << std::endl;
953   cout << "Hexotic execution..." << std::endl;
954   cout << _name << " parameters :" << std::endl;
955   cout << "    " << _name << " Verbosity = " << _hexoticVerbosity << std::endl;
956   cout << "    " << _name << " Max Memory = " << _hexoticMaxMemory << std::endl;
957   cout << "    " << _name << " Segments Min Level = " << _hexesMinLevel << std::endl;
958   cout << "    " << _name << " Segments Max Level = " << _hexesMaxLevel << std::endl;
959   cout << "    " << _name << " Segments Min Size = " << _hexesMinSize << std::endl;
960   cout << "    " << _name << " Segments Max Size = " << _hexesMaxSize << std::endl;
961   cout << "    " << "Hexotic can ignore ridges : " << (_hexoticIgnoreRidges ? "yes":"no") << std::endl;
962   cout << "    " << "Hexotic authorize invalide elements : " << ( _hexoticInvalidElements ? "yes":"no") << std::endl;
963   cout << "    " << _name << " Sharp angle threshold = " << _hexoticSharpAngleThreshold << " degrees" << std::endl;
964   cout << "    " << _name << " Number of threads = " << _hexoticNbProc << std::endl;
965   cout << "    " << _name << " Working directory = \"" << _hexoticWorkingDirectory << "\"" << std::endl;
966   cout << "    " << _name << " Sub. Dom mode = " << _hexoticSdMode << std::endl;
967
968   TCollection_AsciiString run_Hexotic( "mg-hexa.exe" );
969
970   TCollection_AsciiString minl = " --min_level ", maxl = " --max_level ", angle = " --ridge_angle ";
971   TCollection_AsciiString mins = " --min_size ", maxs = " --max_size ";
972   TCollection_AsciiString in   = " --in ",   out  = " --out ";
973   TCollection_AsciiString sizeMap = " --read_sizemap ";
974   TCollection_AsciiString ignoreRidges = " --compute_ridges no ", invalideElements = " --allow_invalid_elements yes ";
975   TCollection_AsciiString subdom = " --components ";
976 #ifndef WIN32
977   TCollection_AsciiString proc = " --max_number_of_threads ";
978 #endif
979   TCollection_AsciiString verb = " --verbose ";
980   TCollection_AsciiString maxmem = " --max_memory ";
981
982   TCollection_AsciiString minLevel, maxLevel, minSize, maxSize, sharpAngle, mode, nbproc, verbosity, maxMemory;
983   minLevel = _hexesMinLevel;
984   maxLevel = _hexesMaxLevel;
985   minSize = _hexesMinSize;
986   maxSize = _hexesMaxSize;
987   sharpAngle = _hexoticSharpAngleThreshold;
988   // Mode translation for mg-tetra 1.1
989   switch ( _hexoticSdMode )
990   {
991     case 1:
992       mode = "outside_skin_only";
993       break;
994     case 2:
995       mode = "outside_components";
996       break;
997     case 3:
998       mode = "all";
999       break;
1000     case 4:
1001       mode = "all --manifold_geometry no";
1002       break;
1003   }
1004   nbproc = _hexoticNbProc;
1005   verbosity = _hexoticVerbosity;
1006   maxMemory = _hexoticMaxMemory;
1007
1008   if (_hexoticIgnoreRidges)
1009     run_Hexotic +=  ignoreRidges;
1010
1011   if (_hexoticInvalidElements)
1012     run_Hexotic +=  invalideElements;
1013
1014   if (_hexesMinSize > 0)
1015     run_Hexotic +=  mins + minSize;
1016
1017   if (_hexesMaxSize > 0)
1018     run_Hexotic +=  maxs + maxSize;
1019
1020   if (_hexesMinLevel > 0)
1021     run_Hexotic +=  minl + minLevel;
1022
1023   if (_hexesMaxLevel > 0)
1024     run_Hexotic +=  maxl + maxLevel;
1025
1026   if (_hexoticSharpAngleThreshold > 0)
1027     run_Hexotic +=  angle + sharpAngle;
1028   
1029   if (_sizeMaps.begin() != _sizeMaps.end())
1030     run_Hexotic += sizeMap + Hexotic_SizeMap_Prefix;
1031
1032   run_Hexotic += in + Hexotic_In + out + Hexotic_Out;
1033   run_Hexotic += subdom + mode;
1034 #ifndef WIN32
1035   run_Hexotic += proc + nbproc;
1036 #endif
1037   run_Hexotic += verb + verbosity;
1038   run_Hexotic += maxmem + maxMemory;
1039
1040   return run_Hexotic.ToCString();
1041 }
1042
1043 // TODO : this is a duplication of some code found in BLSURFPlugin_BLSURF find a proper
1044 // way to share it
1045 TopoDS_Shape HexoticPlugin_Hexotic::entryToShape(std::string entry)
1046 {
1047   MESSAGE("HexoticPlugin_Hexotic::entryToShape "<<entry );
1048   GEOM::GEOM_Object_var aGeomObj;
1049   TopoDS_Shape S = TopoDS_Shape();
1050   SALOMEDS::SObject_var aSObj = myStudy->FindObjectID( entry.c_str() );
1051   if (!aSObj->_is_nil()) {
1052     CORBA::Object_var obj = aSObj->GetObject();
1053     aGeomObj = GEOM::GEOM_Object::_narrow(obj);
1054     aSObj->UnRegister();
1055   }
1056   if ( !aGeomObj->_is_nil() )
1057     S = smeshGen_i->GeomObjectToShape( aGeomObj.in() );
1058   return S;
1059 }
1060
1061 //================================================================================
1062 /*!
1063  * \brief Produces a .mesh file with the size maps informations to give to Hexotic
1064  */
1065 //================================================================================
1066 std::vector<std::string> HexoticPlugin_Hexotic::writeSizeMapFile( std::string sizeMapPrefix )
1067 {
1068   HexoticPlugin_Hypothesis::THexoticSizeMaps::iterator it ;
1069   
1070   // The GMF driver will be used to write the size map file
1071   DriverGMF_Write aWriter;
1072   aWriter.SetSizeMapPrefix( sizeMapPrefix ); 
1073   
1074   std::vector<Control_Pnt> points;
1075   // Iterate on the size maps
1076   for (it=_sizeMaps.begin(); it!=_sizeMaps.end(); it++)
1077   {
1078     // Step 1 : Get the GEOM object entry and the size 
1079     // from the _sizeMaps infos
1080     std::string anEntry = it->first;
1081     double aLocalSize = it->second;
1082     TopoDS_Shape aShape = entryToShape( anEntry );
1083     
1084     // Step 2 : Create the points
1085     createControlPoints( aShape, aLocalSize, points );
1086   }
1087   // Write the .mesh size map file
1088   aWriter.PerformSizeMap(points); 
1089   return aWriter.GetSizeMapFiles();
1090 }
1091
1092 //================================================================================
1093 /*!
1094  * \brief Fills a vector of points from which a size map input file can be written
1095  */
1096 //================================================================================
1097 void HexoticPlugin_Hexotic::createControlPoints( const TopoDS_Shape& aShape, 
1098                                                  const double& theSize, 
1099                                                  std::vector<Control_Pnt>& thePoints )
1100
1101   if ( aShape.ShapeType() == TopAbs_VERTEX )
1102   {
1103     gp_Pnt aPnt = BRep_Tool::Pnt( TopoDS::Vertex(aShape) );
1104     Control_Pnt aControl_Pnt( aPnt, theSize );
1105     thePoints.push_back( aControl_Pnt );
1106   }
1107   if ( aShape.ShapeType() == TopAbs_EDGE )
1108   {
1109     createPointsSampleFromEdge( aShape, theSize, thePoints );  
1110   }
1111   else if ( aShape.ShapeType() == TopAbs_WIRE )
1112   {
1113     TopExp_Explorer Ex;
1114     for (Ex.Init(aShape,TopAbs_EDGE); Ex.More(); Ex.Next()) 
1115     {
1116       createPointsSampleFromEdge( Ex.Current(), theSize, thePoints );
1117     } 
1118   }
1119   else if ( aShape.ShapeType() ==  TopAbs_FACE )
1120   {
1121     createPointsSampleFromFace( aShape, theSize, thePoints ); 
1122   }
1123   else if ( aShape.ShapeType() ==  TopAbs_SOLID )
1124   {
1125     createPointsSampleFromSolid( aShape, theSize, thePoints ); 
1126   }
1127   else if ( aShape.ShapeType() == TopAbs_COMPOUND )
1128   {
1129     TopoDS_Iterator it( aShape );
1130     for(; it.More(); it.Next())
1131     {
1132       createControlPoints( it.Value(), theSize, thePoints );
1133     }
1134   }
1135 }
1136
1137 //================================================================================
1138 /*!
1139  * \brief Fills a vector of points with point samples approximately 
1140  * \brief spaced with a given size
1141  */
1142 //================================================================================
1143 void HexoticPlugin_Hexotic::createPointsSampleFromEdge( const TopoDS_Shape& aShape, 
1144                                                         const double& theSize, 
1145                                                         std::vector<Control_Pnt>& thePoints )
1146 {
1147   double step = theSize;
1148   double first, last;  
1149   Handle( Geom_Curve ) aCurve = BRep_Tool::Curve( TopoDS::Edge( aShape ), first, last );
1150   GeomAdaptor_Curve C ( aCurve );
1151   GCPnts_UniformAbscissa DiscretisationAlgo(C, step , first, last, Precision::Confusion());
1152   int nbPoints = DiscretisationAlgo.NbPoints();
1153   
1154   for ( int i = 1; i <= nbPoints; i++ )
1155   {
1156     double param = DiscretisationAlgo.Parameter( i );
1157     Control_Pnt aPnt;
1158     aCurve->D0( param, aPnt );
1159     aPnt.SetSize(theSize);
1160     thePoints.push_back( aPnt );
1161   }  
1162 }
1163
1164 //================================================================================
1165 /*!
1166  * \brief Fills a vector of points with point samples approximately 
1167  * \brief spaced with a given size
1168  */
1169 //================================================================================
1170 void HexoticPlugin_Hexotic::createPointsSampleFromFace( const TopoDS_Shape& aShape, 
1171                                                         const double& theSize, 
1172                                                         std::vector<Control_Pnt>& thePoints )
1173 {
1174   BRepMesh_IncrementalMesh M(aShape, 0.01, Standard_True);
1175   TopLoc_Location aLocation;
1176   TopoDS_Face aFace = TopoDS::Face(aShape);
1177
1178   // Triangulate the face
1179   Handle(Poly_Triangulation) aTri = BRep_Tool::Triangulation (aFace, aLocation);
1180   
1181   // Get the transformation associated to the face location
1182   gp_Trsf aTrsf = aLocation.Transformation();
1183   
1184   // Get triangles
1185   int nbTriangles = aTri->NbTriangles();
1186   Poly_Array1OfTriangle triangles(1,nbTriangles);
1187   triangles=aTri->Triangles();
1188   
1189   // GetNodes
1190   int nbNodes = aTri->NbNodes();
1191   TColgp_Array1OfPnt nodes(1,nbNodes);
1192   nodes = aTri->Nodes();
1193
1194   // Iterate on triangles and subdivide them
1195   for(int i=1; i<=nbTriangles; i++)
1196   {
1197      Poly_Triangle aTriangle = triangles.Value(i);
1198      gp_Pnt p1 = nodes.Value(aTriangle.Value(1));
1199      gp_Pnt p2 = nodes.Value(aTriangle.Value(2));
1200      gp_Pnt p3 = nodes.Value(aTriangle.Value(3));
1201      
1202      p1.Transform(aTrsf);
1203      p2.Transform(aTrsf);
1204      p3.Transform(aTrsf);
1205      
1206      subdivideTriangle(p1, p2, p3, theSize, thePoints);  
1207   }
1208 }
1209
1210 //================================================================================
1211 /*!
1212  * \brief Fills a vector of points with point samples approximately 
1213  * \brief spaced with a given size
1214  */
1215 //================================================================================
1216 void HexoticPlugin_Hexotic::createPointsSampleFromSolid( const TopoDS_Shape& aShape, 
1217                                                          const double& theSize, 
1218                                                          std::vector<Control_Pnt>& thePoints )
1219 {
1220   // Compute the bounding box
1221   double Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
1222   Bnd_Box B;               
1223   BRepBndLib::Add(aShape, B);
1224   B.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
1225   
1226   // Create the points
1227   double step = theSize;
1228   
1229   for ( double x=Xmin; x-Xmax<Precision::Confusion(); x=x+step )
1230   {
1231     for ( double y=Ymin; y-Ymax<Precision::Confusion(); y=y+step )
1232     {
1233       // Step1 : generate the Zmin -> Zmax line
1234       gp_Pnt startPnt(x, y, Zmin);
1235       gp_Pnt endPnt(x, y, Zmax);
1236       gp_Vec aVec(startPnt, endPnt);
1237       gp_Lin aLine(startPnt, aVec);
1238       double endParam = Zmax - Zmin;
1239       
1240       // Step2 : for each face of aShape:
1241       std::set<double> intersections;
1242       std::set<double>::iterator it = intersections.begin();
1243       
1244       TopExp_Explorer Ex;
1245       for (Ex.Init(aShape,TopAbs_FACE); Ex.More(); Ex.Next()) 
1246       { 
1247         // check if there is an intersection
1248         IntCurvesFace_Intersector anIntersector(TopoDS::Face(Ex.Current()), Precision::Confusion());
1249         anIntersector.Perform(aLine, 0, endParam);
1250         
1251         // get the intersection's parameter and store it
1252         int nbPoints = anIntersector.NbPnt();
1253         for(int i = 0 ; i < nbPoints ; i++ )
1254         {
1255           it = intersections.insert( it, anIntersector.WParameter(i+1) );
1256         }
1257       }
1258       // Step3 : go through the line chunk by chunk 
1259       if ( intersections.begin() != intersections.end() )
1260       {
1261         std::set<double>::iterator intersectionsIterator=intersections.begin();
1262         double first = *intersectionsIterator;
1263         intersectionsIterator++;
1264         bool innerPoints = true; 
1265         for ( ; intersectionsIterator!=intersections.end() ; intersectionsIterator++ )
1266         {
1267           double second = *intersectionsIterator;
1268           if ( innerPoints )
1269           {
1270             // If the last chunk was outside of the shape or this is the first chunk
1271             // add the points in the range [first, second] to the points vector
1272             double localStep = (second -first) / ceil( (second - first) / step );
1273             for ( double z = Zmin + first; z < Zmin + second; z = z + localStep )
1274             {
1275               thePoints.push_back(Control_Pnt( x, y, z, theSize ));
1276             }
1277             thePoints.push_back(Control_Pnt( x, y, Zmin + second, theSize ));
1278           }
1279           first = second;
1280           innerPoints = !innerPoints;
1281         }
1282       }
1283     }
1284   }
1285 }
1286
1287 //================================================================================
1288 /*!
1289  * \brief Subdivides a triangle until it reaches a certain size (recursive function)
1290  */
1291 //================================================================================
1292 void HexoticPlugin_Hexotic::subdivideTriangle( const gp_Pnt& p1, 
1293                                                const gp_Pnt& p2, 
1294                                                const gp_Pnt& p3, 
1295                                                const double& theSize, 
1296                                                std::vector<Control_Pnt>& thePoints)
1297 {
1298   // Size threshold to stop subdividing
1299   // This value ensures that two control points are distant no more than 2*theSize
1300   // as shown below
1301   //
1302   // The greater distance D of the mass center M to each Edge is 1/3 * Median 
1303   // and Median < sqrt(3/4) * a  where a is the greater side (by using Apollonius' thorem). 
1304   // So D < 1/3 * sqrt(3/4) * a and if a < sqrt(3) * S then D < S/2
1305   // and the distance between two mass centers of two neighbouring triangles 
1306   // sharing an edge is < 2 * 1/2 * S = S
1307   // If the traingles share a Vertex and no Edge the distance of the mass centers 
1308   // to the Vertices is 2*D < S so the mass centers are distant of less than 2*S 
1309   
1310   double threshold = sqrt( 3. ) * theSize;
1311   
1312   if ( (p1.Distance(p2) > threshold ||
1313         p2.Distance(p3) > threshold ||
1314         p3.Distance(p1) > threshold))
1315   { 
1316     std::vector<gp_Pnt> midPoints = computePointsForSplitting(p1, p2, p3);
1317  
1318     subdivideTriangle( midPoints[0], midPoints[1], midPoints[2], theSize, thePoints );
1319     subdivideTriangle( midPoints[0], p2, midPoints[1], theSize, thePoints );
1320     subdivideTriangle( midPoints[2], midPoints[1], p3, theSize, thePoints );
1321     subdivideTriangle( p1, midPoints[0], midPoints[2], theSize, thePoints );
1322   }
1323   else
1324   {
1325     double x = (p1.X() + p2.X() + p3.X()) / 3 ;
1326     double y = (p1.Y() + p2.Y() + p3.Y()) / 3 ;
1327     double z = (p1.Z() + p2.Z() + p3.Z()) / 3 ;
1328     
1329     Control_Pnt massCenter( x ,y ,z, theSize );
1330     thePoints.push_back( massCenter );
1331   }
1332 }
1333
1334 //================================================================================
1335 /*!
1336  * \brief Returns the appropriate points for splitting a triangle
1337  * \brief the tangency points of the incircle are used in order to have mostly
1338  * \brief well-shaped sub-triangles
1339  */
1340 //================================================================================
1341 std::vector<gp_Pnt> HexoticPlugin_Hexotic::computePointsForSplitting( const gp_Pnt& p1, 
1342                                                                       const gp_Pnt& p2, 
1343                                                                       const gp_Pnt& p3 )
1344 {
1345   std::vector<gp_Pnt> midPoints;
1346   //Change coordinates
1347   gp_Trsf Trsf_1;            // Identity transformation
1348   gp_Ax3 reference_system(gp::Origin(), gp::DZ(), gp::DX());   // OXY
1349  
1350   gp_Vec Vx(p1, p3);
1351   gp_Vec Vaux(p1, p2);
1352   gp_Dir Dx(Vx);
1353   gp_Dir Daux(Vaux);
1354   gp_Dir Dz = Dx.Crossed(Daux);
1355   gp_Ax3 current_system(p1, Dz, Dx);
1356   
1357   Trsf_1.SetTransformation( reference_system, current_system );
1358   
1359   gp_Pnt A = p1.Transformed(Trsf_1);
1360   gp_Pnt B = p2.Transformed(Trsf_1);
1361   gp_Pnt C = p3.Transformed(Trsf_1);
1362   
1363   double a =  B.Distance(C) ;
1364   double b =  A.Distance(C) ;
1365   double c =  B.Distance(A) ;
1366   
1367   // Incenter coordinates
1368   // see http://mathworld.wolfram.com/Incenter.html
1369   double Xi = ( b*B.X() + c*C.X() ) / ( a + b + c );
1370   double Yi = ( b*B.Y() ) / ( a + b + c );
1371   gp_Pnt Center(Xi, Yi, 0);
1372   
1373   // Calculate the tangency points of the incircle
1374   gp_Pnt T1 = tangencyPoint( A, B, Center);
1375   gp_Pnt T2 = tangencyPoint( B, C, Center);
1376   gp_Pnt T3 = tangencyPoint( C, A, Center);
1377   
1378   gp_Pnt p1_2 = T1.Transformed(Trsf_1.Inverted());
1379   gp_Pnt p2_3 = T2.Transformed(Trsf_1.Inverted());
1380   gp_Pnt p3_1 = T3.Transformed(Trsf_1.Inverted());
1381
1382   midPoints.push_back(p1_2);
1383   midPoints.push_back(p2_3);
1384   midPoints.push_back(p3_1);
1385   
1386   return midPoints;
1387 }
1388
1389 //================================================================================
1390 /*!
1391  * \brief Computes the tangency points of the circle of center Center with
1392  * \brief the straight line (p1 p2)
1393  */
1394 //================================================================================
1395 gp_Pnt HexoticPlugin_Hexotic::tangencyPoint(const gp_Pnt& p1,
1396                                             const gp_Pnt& p2,
1397                                             const gp_Pnt& Center)
1398 {
1399   double Xt = 0;
1400   double Yt = 0;
1401   
1402   // The tangency point is the intersection of the straight line (p1 p2)
1403   // and the straight line (Center T) which is orthogonal to (p1 p2)
1404   if ( fabs(p1.X() - p2.X()) <= Precision::Confusion() )
1405   {
1406     Xt=p1.X();     // T is on (p1 p2)
1407     Yt=Center.Y(); // (Center T) is orthogonal to (p1 p2)
1408   }
1409   else if ( fabs(p1.Y() - p2.Y()) <= Precision::Confusion() )
1410   {
1411     Yt=p1.Y();     // T is on (p1 p2) 
1412     Xt=Center.X(); // (Center T) is orthogonal to (p1 p2)
1413   }
1414   else
1415   {
1416     // First straight line coefficients (equation y=a*x+b)
1417     double a = (p2.Y() - p1.Y()) / (p2.X() - p1.X())  ;
1418     double b = p1.Y() - a*p1.X();         // p1 is on this straight line
1419     
1420     // Second straight line coefficients (equation y=c*x+d)
1421     double c = -1 / a;                    // The 2 lines are orthogonal
1422     double d = Center.Y() - c*Center.X(); // Center is on this straight line
1423     
1424     Xt = (d - b) / (a - c);
1425     Yt = a*Xt + b;
1426   }
1427   
1428   return gp_Pnt( Xt, Yt, 0 );
1429 }
1430
1431 //=============================================================================
1432 /*!
1433  * Here we are going to use the Hexotic mesher
1434  */
1435 //=============================================================================
1436
1437 bool HexoticPlugin_Hexotic::Compute(SMESH_Mesh&          aMesh,
1438                                     const TopoDS_Shape& aShape)
1439 {
1440   _compute_canceled = false;
1441   bool Ok = true;
1442   SMESHDS_Mesh* meshDS = aMesh.GetMeshDS();
1443   TCollection_AsciiString hexahedraMessage;
1444
1445   if (_iShape == 0 && _nbShape == 0) {
1446     _nbShape = countShape( meshDS, TopAbs_SOLID );  // we count the number of shapes
1447   }
1448
1449   // to prevent from displaying error message after computing,
1450   // SetIsAlwaysComputed( true ) to empty sub-meshes
1451   vector< SMESH_subMesh* > subMeshesAlwaysComp;
1452   for ( int i = 0; i < _nbShape; ++i )
1453     if ( SMESH_subMesh* sm = aMesh.GetSubMeshContaining( aShape ))
1454     {
1455       SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(/*includeSelf=*/true,
1456                                                                /*complexShapeFirst=*/false);
1457       while ( smIt->more() )
1458       {
1459         sm = smIt->next();
1460         if ( !sm->IsMeshComputed() )
1461         {
1462           sm->SetIsAlwaysComputed( true );
1463           subMeshesAlwaysComp.push_back( sm );
1464         }
1465       }
1466     }
1467
1468   _iShape++;
1469
1470   if (_iShape == _nbShape ) {
1471
1472     // create bounding box for each shape of the compound
1473
1474     int iShape = 0;
1475     TopoDS_Shape *tabShape;
1476     double **tabBox;
1477
1478     tabShape = new TopoDS_Shape[_nbShape];
1479     tabBox   = new double*[_nbShape];
1480     for (int i=0; i<_nbShape; i++)
1481       tabBox[i] = new double[6];
1482     double Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
1483
1484     TopExp_Explorer expBox (meshDS->ShapeToMesh(), TopAbs_SOLID);
1485     for (; expBox.More(); expBox.Next()) {
1486       tabShape[iShape] = expBox.Current();
1487       Bnd_Box BoundingBox;
1488       BRepBndLib::Add(expBox.Current(), BoundingBox);
1489       BoundingBox.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
1490       tabBox[iShape][0] = Xmin; tabBox[iShape][1] = Xmax;
1491       tabBox[iShape][2] = Ymin; tabBox[iShape][3] = Ymax;
1492       tabBox[iShape][4] = Zmin; tabBox[iShape][5] = Zmax;
1493       iShape++;
1494     }
1495
1496     SetParameters(_hypothesis);
1497
1498 //     TCollection_AsciiString aTmpDir = getTmpDir();
1499     TCollection_AsciiString aTmpDir = _hexoticWorkingDirectory.c_str();
1500 #ifdef WIN32
1501     if ( aTmpDir.Value(aTmpDir.Length()) != '\\' ) aTmpDir += '\\';
1502 #else
1503     if ( aTmpDir.Value(aTmpDir.Length()) != '/' ) aTmpDir += '/';
1504 #endif
1505     TCollection_AsciiString Hexotic_In(""), Hexotic_Out, Hexotic_SizeMap_Prefix;
1506     TCollection_AsciiString modeFile_In( "chmod 666 " ), modeFile_Out( "chmod 666 " );
1507     TCollection_AsciiString aLogFileName = aTmpDir + "Hexotic"+getSuffix()+".log";    // log
1508
1509     std::map <int,int> aSmdsToHexoticIdMap;
1510     std::map <int,const SMDS_MeshNode*> aHexoticIdToNodeMap;
1511
1512     Hexotic_Out = aTmpDir + "Hexotic"+getSuffix()+"_Out.mesh";
1513 #ifdef WITH_BLSURFPLUGIN
1514     bool defaultInputFile = true;
1515     if (_blsurfHypo && !_blsurfHypo->GetQuadAllowed()) {
1516       Hexotic_In = TCollection_AsciiString(_blsurfHypo->GetGMFFile().c_str());
1517       if (Hexotic_In != "")
1518         defaultInputFile = false;
1519     }
1520     if (defaultInputFile) {
1521 #endif
1522       Hexotic_In  = aTmpDir + "Hexotic"+getSuffix()+"_In.mesh";
1523       removeHexoticFiles(Hexotic_In, Hexotic_Out);
1524       splitQuads(aMesh); // quadrangles are no longer acceptable as input
1525       cout << std::endl;
1526       cout << "Creating Hexotic input mesh file : " << Hexotic_In << std::endl;
1527       aMesh.ExportGMF(Hexotic_In.ToCString(), meshDS, true);
1528 #ifdef WITH_BLSURFPLUGIN
1529     }
1530     else {
1531       removeFile( Hexotic_Out );
1532     }
1533 #endif
1534     
1535     Hexotic_SizeMap_Prefix = aTmpDir + "Hexotic_SizeMap" + getSuffix();
1536     std::vector<std::string> sizeMapFiles = writeSizeMapFile( Hexotic_SizeMap_Prefix.ToCString() );
1537     
1538     std::string run_Hexotic = getHexoticCommand(Hexotic_In, Hexotic_Out, Hexotic_SizeMap_Prefix);
1539     run_Hexotic += std::string(" 1> ") + aLogFileName.ToCString();  // dump into file
1540
1541     cout << std::endl;
1542     cout << "Hexotic command : " << run_Hexotic << std::endl;
1543
1544 #ifndef WIN32    
1545     modeFile_In += Hexotic_In;
1546     system( modeFile_In.ToCString() );
1547 #endif
1548     aSmdsToHexoticIdMap.clear();
1549     aHexoticIdToNodeMap.clear();
1550
1551     MESSAGE("HexoticPlugin_Hexotic::Compute");
1552
1553     int status = system( run_Hexotic.data() );
1554
1555     // --------------
1556     // read a result
1557     // --------------
1558
1559     std::ifstream fileRes( Hexotic_Out.ToCString() );
1560     modeFile_Out += Hexotic_Out;
1561     system( modeFile_Out.ToCString() );
1562     if ( ! fileRes.fail() ) {
1563       Ok = readResult( Hexotic_Out.ToCString(),
1564                        this,
1565                        meshDS, _nbShape, tabShape, tabBox );
1566       if(Ok) {
1567 /*********************
1568 // TODO: Detect and remove elements in holes in case of sd mode = 4
1569       // Remove previous nodes and elements
1570       SMDS_ElemIteratorPtr itElement = meshDS->elementsIterator();
1571       SMDS_NodeIteratorPtr itNode = meshDS->nodesIterator();
1572     
1573       while ( itElement->more() )
1574         meshDS->RemoveElement( itElement->next() );
1575       while ( itNode->more() )
1576         meshDS->RemoveNode( itNode->next() );
1577   
1578       SMESH_ComputeErrorPtr myError = aMesh.GMFToMesh(Hexotic_Out.ToCString());
1579       if (myError)
1580 */
1581         hexahedraMessage = "success";
1582         #ifndef _DEBUG_
1583         removeFile(Hexotic_Out);
1584         removeFile(Hexotic_In);
1585         removeFile(aLogFileName);
1586         for( int i=0; i<sizeMapFiles.size(); i++)
1587         {
1588           removeFile( TCollection_AsciiString( sizeMapFiles[i].c_str() ) );
1589         }
1590         #endif
1591       }
1592       else {
1593         hexahedraMessage = "failed"; 
1594       }
1595     }
1596     else {
1597       hexahedraMessage = "failed";
1598       cout << "Problem with Hexotic output file " << Hexotic_Out.ToCString() << std::endl;
1599       Ok = false;
1600       // analyse log file
1601       SMESH_File logFile( aLogFileName.ToCString() );
1602       if ( !logFile.eof() )
1603       {
1604         char msgLic[] = " Dlim ";
1605         const char* fileBeg = logFile.getPos(), *fileEnd = fileBeg + logFile.size();
1606         if ( std::search( fileBeg, fileEnd, msgLic, msgLic+strlen(msgLic)) != fileEnd )
1607           error("Licence problems.");
1608       }
1609 #ifndef WIN32
1610       if ( status > 0 && WEXITSTATUS(status) == 127 )
1611         error("hexotic: command not found");
1612 #else
1613       int err = errno;
1614       if ( status == 0 && err == ENOENT ) {
1615         error("hexotic: command not found");
1616       }
1617 #endif
1618     }
1619     cout << "Hexahedra meshing " << hexahedraMessage << std::endl;
1620     cout << std::endl;
1621
1622     // restore "always computed" flag of sub-meshes (0022127)
1623     for  ( size_t iSM = 0; iSM < subMeshesAlwaysComp.size(); ++iSM )
1624       subMeshesAlwaysComp[ iSM ]->SetIsAlwaysComputed( false );
1625
1626     delete [] tabShape;
1627     for (int i=0; i<_nbShape; i++)
1628       delete [] tabBox[i];
1629     delete [] tabBox;
1630     _nbShape = 0;
1631     _iShape  = 0;
1632   }
1633   if(_compute_canceled)
1634     return error(SMESH_Comment("interruption initiated by user"));
1635   return Ok;
1636 }
1637
1638 //=============================================================================
1639 /*!
1640  * \brief Computes mesh without geometry
1641  *  \param aMesh - the mesh
1642  *  \param aHelper - helper that must be used for adding elements to \aaMesh
1643  *  \retval bool - is a success
1644  *
1645  * The method is called if ( !aMesh->HasShapeToMesh() )
1646  */
1647 //=============================================================================
1648
1649 bool HexoticPlugin_Hexotic::Compute(SMESH_Mesh & aMesh, SMESH_MesherHelper* aHelper)
1650 {
1651   _compute_canceled = false;
1652 /*
1653   SMESH_ComputeErrorPtr myError = SMESH_ComputeError::New();
1654 */
1655   bool Ok = true;
1656   TCollection_AsciiString hexahedraMessage;
1657
1658   SetParameters(_hypothesis);
1659
1660   TCollection_AsciiString aTmpDir = _hexoticWorkingDirectory.c_str();//getTmpDir();
1661   TCollection_AsciiString Hexotic_In, Hexotic_Out, Hexotic_SizeMap_Prefix;
1662   TCollection_AsciiString modeFile_In( "chmod 666 " ), modeFile_Out( "chmod 666 " );
1663   TCollection_AsciiString aLogFileName = aTmpDir + "Hexotic"+getSuffix()+".log";    // log
1664
1665   std::map <int,int> aSmdsToHexoticIdMap;
1666   std::map <int,const SMDS_MeshNode*> aHexoticIdToNodeMap;
1667
1668   Hexotic_In  = aTmpDir + "Hexotic"+getSuffix()+"_In.mesh";
1669   Hexotic_Out = aTmpDir + "Hexotic"+getSuffix()+"_Out.mesh";
1670   Hexotic_SizeMap_Prefix = aTmpDir + "Hexotic_SizeMap";
1671  
1672   
1673   std::vector<std::string> sizeMapFiles = writeSizeMapFile( Hexotic_SizeMap_Prefix.ToCString() );
1674
1675   std::string run_Hexotic = getHexoticCommand(Hexotic_In, Hexotic_Out, Hexotic_SizeMap_Prefix);
1676   run_Hexotic += std::string(" 1> ") + aLogFileName.ToCString();  // dump into file
1677
1678   removeHexoticFiles(Hexotic_In, Hexotic_Out);
1679
1680   splitQuads(aMesh); // quadrangles are no longer acceptable as input
1681
1682   cout << std::endl;
1683   cout << "Creating Hexotic input mesh file : " << Hexotic_In << std::endl;
1684   aMesh.ExportGMF(Hexotic_In.ToCString(), aHelper->GetMeshDS());
1685 #ifndef WIN32    
1686   modeFile_In += Hexotic_In;
1687   system( modeFile_In.ToCString() );
1688 #endif
1689   aSmdsToHexoticIdMap.clear();
1690   aHexoticIdToNodeMap.clear();
1691
1692   MESSAGE("HexoticPlugin_Hexotic::Compute");
1693
1694   cout << std::endl;
1695   cout << "Hexotic command : " << run_Hexotic << std::endl;
1696   system( run_Hexotic.data() );
1697
1698   // --------------
1699   // read a result
1700   // --------------
1701
1702   std::ifstream fileRes( Hexotic_Out.ToCString() );
1703   modeFile_Out += Hexotic_Out;
1704   system( modeFile_Out.ToCString() );
1705   if ( ! fileRes.fail() ) {
1706     Ok = readResult( Hexotic_Out.ToCString(),
1707                      this,
1708                      aHelper );
1709     if(Ok)
1710 /*
1711     // Remove previous nodes and elements
1712     SMDS_ElemIteratorPtr itElement = aHelper->GetMeshDS()->elementsIterator();
1713     SMDS_NodeIteratorPtr itNode = aHelper->GetMeshDS()->nodesIterator();
1714     
1715     while ( itElement->more() )
1716       aHelper->GetMeshDS()->RemoveElement( itElement->next() );
1717     while ( itNode->more() )
1718       aHelper->GetMeshDS()->RemoveNode( itNode->next() );
1719
1720     // Import GMF mesh
1721     myError = aMesh.GMFToMesh(Hexotic_Out.ToCString());
1722     
1723     itElement = aHelper->GetMeshDS()->elementsIterator();
1724     itNode = aHelper->GetMeshDS()->nodesIterator();
1725
1726     // Assign nodes and elements to the pseudo shape
1727     while ( itNode->more() )
1728       aHelper->GetMeshDS()->SetNodeInVolume(itNode->next(), 1);
1729     while ( itElement->more() )
1730       aHelper->GetMeshDS()->SetMeshElementOnShape(itElement->next(), 1);
1731
1732     if(myError->IsOK())
1733 */
1734       hexahedraMessage = "success";
1735     else
1736       hexahedraMessage = "failed";
1737   }
1738   else {
1739 /*
1740     myError->myName = COMPERR_EXCEPTION;
1741 */
1742     hexahedraMessage = "failed";
1743     cout << "Problem with Hexotic output file " << Hexotic_Out << std::endl;
1744     // analyse log file
1745     SMESH_File logFile( aLogFileName.ToCString() );
1746     if ( !logFile.eof() )
1747     {
1748       char msgLic[] = " Dlim ";
1749       const char* fileBeg = logFile.getPos(), *fileEnd = fileBeg + logFile.size();
1750       if ( std::search( fileBeg, fileEnd, msgLic, msgLic+strlen(msgLic)) != fileEnd )
1751         return error("Licence problems.");
1752     }
1753     return error(SMESH_Comment("Problem with Hexotic output file ")<<Hexotic_Out);
1754   }
1755   cout << "Hexahedra meshing " << hexahedraMessage << std::endl;
1756   cout << std::endl;
1757
1758   if(_compute_canceled)
1759     return error(SMESH_Comment("interruption initiated by user"));
1760   removeFile(Hexotic_Out);
1761   removeFile(Hexotic_In);
1762   removeFile(aLogFileName);
1763   for( int i=0; i<sizeMapFiles.size(); i++)
1764   {
1765     removeFile( TCollection_AsciiString(sizeMapFiles[i].c_str()) );
1766   }
1767   return Ok;
1768 /*
1769   return myError->IsOK();
1770 */
1771 }
1772
1773 //=============================================================================
1774 /*!
1775  *
1776  */
1777 //=============================================================================
1778
1779 bool HexoticPlugin_Hexotic::Evaluate(SMESH_Mesh& aMesh,
1780                                      const TopoDS_Shape& aShape,
1781                                      MapShapeNbElems& aResMap)
1782 {
1783   std::vector<int> aResVec(SMDSEntity_Last);
1784   for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
1785   SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
1786   aResMap.insert(std::make_pair(sm,aResVec));
1787
1788   SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1789   smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,"Evaluation is not implemented",this));
1790
1791   return true;
1792 }
1793
1794 void HexoticPlugin_Hexotic::CancelCompute()
1795 {
1796   _compute_canceled = true;
1797 #ifdef WIN32
1798 #else
1799   TCollection_AsciiString aTmpDir = _hexoticWorkingDirectory.c_str(); //getTmpDir();
1800   TCollection_AsciiString Hexotic_In = aTmpDir + "Hexotic_In.mesh";
1801   TCollection_AsciiString cmd = TCollection_AsciiString("ps ux | grep ") + Hexotic_In;
1802   cmd += TCollection_AsciiString(" | grep -v grep | awk '{print $2}' | xargs kill -9 > /dev/null 2>&1");
1803   system( cmd.ToCString() );
1804 #endif
1805 }