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