Salome HOME
0052232: TC7.3.0: unexpected result in the mesh computation
[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 #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 : readResult
374 //purpose  : Read GMF file in case of a mesh with geometry
375 //=======================================================================
376
377 static bool readResult(std::string         theFile,
378                        HexoticPlugin_Hexotic*  theAlgo,
379                        SMESHDS_Mesh*       theMesh,
380                        const int           nbShape,
381                        const TopoDS_Shape* tabShape,
382                        double**            tabBox)
383 {
384   // ---------------------------------
385   // Optimisation of the plugin ...
386   // Retrieve the correspondance edge --> shape
387   // (which is very costly) only when user
388   // has defined at least one group of edges
389   // which should be rare for a 3d mesh !
390   // ---------------------------------
391   
392   bool retrieve_edges = false;
393   const std::set<SMESHDS_GroupBase*>& aGroups = theMesh->GetGroups();
394   set<SMESHDS_GroupBase*>::const_iterator GrIt = aGroups.begin();
395   for (; GrIt != aGroups.end(); GrIt++)
396     {
397       SMESHDS_GroupBase* aGrp = *GrIt;
398       if ( !aGrp )
399         continue;
400       if ( aGrp->GetType() == SMDSAbs_Edge )
401         {
402           retrieve_edges = true;
403           break;
404         }
405     }
406   
407   // ---------------------------------
408   // Read generated elements and nodes
409   // ---------------------------------
410
411   TopoDS_Shape aShape;
412   TopoDS_Vertex aVertex;
413   std::string token;
414   int EndOfFile = 0, nbElem = 0, nField = 9, nbRef = 0;
415   int aHexoticNodeID = 0, shapeID, hexoticShapeID;
416   const int IdShapeRef = 2;
417   int *tabID, *tabRef, *nodeAssigne;
418   bool *tabDummy, hasDummy = false;
419   double epsilon = Precision::Confusion();
420   std::map <std::string,int> mapField;
421   SMDS_MeshNode** HexoticNode;
422   TopoDS_Shape *tabCorner, *tabEdge;
423
424   const int nbDomains = countShape( theMesh, TopAbs_SHELL );
425   const int holeID = -1;
426
427   // tabID    = new int[nbShape];
428   tabID    = new int[nbDomains];
429   tabRef   = new int[nField];
430   tabDummy = new bool[nField];
431
432   for (int i=0; i<nbDomains; i++)
433     tabID[i] = 0;
434   if ( nbDomains == 1 )
435     tabID[0] = theMesh->ShapeToIndex( tabShape[0] );
436
437   mapField["MeshVersionFormatted"] = 0; tabRef[0] = 0; tabDummy[0] = false;
438   mapField["Dimension"]            = 1; tabRef[1] = 0; tabDummy[1] = false;
439   mapField["Vertices"]             = 2; tabRef[2] = 3; tabDummy[2] = true;
440   mapField["Corners"]              = 3; tabRef[3] = 1; tabDummy[3] = false;
441   mapField["Edges"]                = 4; tabRef[4] = 2; tabDummy[4] = true;
442   mapField["Ridges"]               = 5; tabRef[5] = 1; tabDummy[5] = false;
443   mapField["Quadrilaterals"]       = 6; tabRef[6] = 4; tabDummy[6] = true;
444   mapField["Hexahedra"]            = 7; tabRef[7] = 8; tabDummy[7] = true;
445   mapField["End"]                  = 8; tabRef[8] = 0; tabDummy[0] = false;
446
447   SMDS_NodeIteratorPtr itOnHexoticInputNode = theMesh->nodesIterator();
448   while ( itOnHexoticInputNode->more() )
449     theMesh->RemoveNode( itOnHexoticInputNode->next() );
450
451   int nbVertices   = getNbShape(theFile, "Vertices");
452   int nbCorners    = getNbShape(theFile, "Corners", countShape( theMesh, TopAbs_VERTEX ));
453   int nbShapeEdge  = countShape( theMesh, TopAbs_EDGE );
454
455   tabCorner   = new TopoDS_Shape[ nbCorners ];
456   tabEdge     = new TopoDS_Shape[ nbShapeEdge ];
457   nodeAssigne = new int[ nbVertices + 1 ];
458   HexoticNode = new SMDS_MeshNode*[ nbVertices + 1 ];
459
460   getShape(theMesh, TopAbs_VERTEX, tabCorner);
461   getShape(theMesh, TopAbs_EDGE,   tabEdge);
462
463   MESSAGE("Read " << theFile << " file");
464   std::ifstream fileRes(theFile.c_str());
465   ASSERT(fileRes);
466
467   while ( EndOfFile == 0  ) {
468     int dummy;
469     fileRes >> token;
470
471     if (mapField.count(token)) {
472       nField   = mapField[token];
473       nbRef    = tabRef[nField];
474       hasDummy = tabDummy[nField];
475     }
476     else {
477       nField = -1;
478       nbRef = 0;
479     }
480
481     nbElem = 0;
482     if ( nField < (mapField.size() - 1) && nField >= 0 )
483       fileRes >> nbElem;
484
485     switch (nField) {
486       case 0: { // "MeshVersionFormatted"
487         MESSAGE(token << " " << nbElem);
488         break;
489       }
490       case 1: { // "Dimension"
491         MESSAGE("Mesh dimension " << nbElem << "D");
492         break;
493       }
494       case 2: { // "Vertices"
495         MESSAGE("Read " << nbElem << " " << token);
496         int aHexoticID;
497         double *coord;
498         SMDS_MeshNode * aHexoticNode;
499
500         coord = new double[nbRef];
501         for ( int iElem = 0; iElem < nbElem; iElem++ ) {
502           if(theAlgo->computeCanceled())
503             {
504               return false;
505             }
506           aHexoticID = iElem + 1;
507           for ( int iCoord = 0; iCoord < 3; iCoord++ )
508             fileRes >> coord[ iCoord ];
509           fileRes >> dummy;
510           aHexoticNode = theMesh->AddNode(coord[0], coord[1], coord[2]);
511           HexoticNode[ aHexoticID ] = aHexoticNode;
512           nodeAssigne[ aHexoticID ] = 0;
513         }
514         delete [] coord;
515         break;
516       }
517       case 3: // "Corners"
518       case 4: // "Edges"
519       case 5: // "Ridges"
520       case 6: // "Quadrilaterals"
521       case 7: { // "Hexahedra"
522         MESSAGE("Read " << nbElem << " " << token);
523         SMDS_MeshNode** node;
524         int nodeDim, *nodeID;
525         SMDS_MeshElement * aHexoticElement = 0;
526
527         node   = new SMDS_MeshNode*[ nbRef ];
528         nodeID = new int[ nbRef ];
529         for ( int iElem = 0; iElem < nbElem; iElem++ ) {
530           if(theAlgo->computeCanceled())
531             {
532               return false;
533             }
534           for ( int iRef = 0; iRef < nbRef; iRef++ ) {
535             fileRes >> aHexoticNodeID;                          // read nbRef aHexoticNodeID
536             node[ iRef ]   = HexoticNode[ aHexoticNodeID ];
537             nodeID[ iRef ] = aHexoticNodeID;
538           }
539           if ( hasDummy )
540             fileRes >> dummy;
541           switch (nField) {
542             case 3: { // "Corners"
543               nodeDim = 1;
544               gp_Pnt HexoticPnt ( node[0]->X(), node[0]->Y(), node[0]->Z() );
545               for ( int i=0; i<nbElem; i++ ) {
546                 aVertex = TopoDS::Vertex( tabCorner[i] );
547                 gp_Pnt aPnt = BRep_Tool::Pnt( aVertex );
548                 if ( aPnt.Distance( HexoticPnt ) < epsilon )
549                   break;
550               }
551               break;
552             }
553             case 4: { // "Edges"
554               nodeDim = 2;
555               aHexoticElement = theMesh->AddEdge( node[0], node[1] );
556               int iNode = 1;
557               if ( nodeAssigne[ nodeID[0] ] == 0 || nodeAssigne[ nodeID[0] ] == 2 )
558                 iNode = 0;
559               if(retrieve_edges)
560                 shapeID = findEdge( node[iNode], theMesh, nbShapeEdge, tabEdge );
561               else
562                 shapeID = 0;
563               break;
564             }
565             case 5: { // "Ridges"
566               break;
567             }
568             case 6: { // "Quadrilaterals"
569               nodeDim = 3;
570               aHexoticElement = theMesh->AddFace( node[0], node[1], node[2], node[3] );
571               shapeID = dummy;
572               break;
573             }
574             case 7: { // "Hexahedra"
575               nodeDim = 4;
576               if ( nbDomains > 1 ) {
577                 hexoticShapeID = dummy - IdShapeRef;
578                 if ( tabID[ hexoticShapeID ] == 0 ) {
579                   aShape = findShape(node, aShape, tabShape, tabBox, nbShape);
580                   shapeID = aShape.IsNull() ? holeID : theMesh->ShapeToIndex( aShape );
581                   tabID[ hexoticShapeID ] = shapeID;
582                 }
583                 else
584                   shapeID = tabID[ hexoticShapeID ];
585                 if ( iElem == (nbElem - 1) ) {
586                   int shapeAssociated = 0;
587                   for ( int i=0; i<nbDomains; i++ ) {
588                     if (tabID[i] > 0 )
589                       shapeAssociated += 1;
590                   }
591                   if ( shapeAssociated != nbShape )
592                     printWarning(nbShape, "domains", shapeAssociated);
593                 }
594               }
595               else {
596                 shapeID = tabID[0];
597               }
598               if ( shapeID != holeID )
599                 aHexoticElement = theMesh->AddVolume( node[0], node[3], node[2], node[1], node[4], node[7], node[6], node[5] );
600               break;
601             }
602           } // switch (nField)
603
604           if ( token != "Ridges" && ( shapeID > 0 || token == "Corners")) {
605             for ( int i=0; i<nbRef; i++ ) {
606               if ( nodeAssigne[ nodeID[i] ] == 0 ) {
607                 if      ( token == "Corners" )        theMesh->SetNodeOnVertex( node[0], aVertex );
608                 else if ( token == "Edges" )          theMesh->SetNodeOnEdge( node[i], shapeID );
609                 else if ( token == "Quadrilaterals" ) theMesh->SetNodeOnFace( node[i], shapeID );
610                 else if ( token == "Hexahedra" )      theMesh->SetNodeInVolume( node[i], shapeID );
611                 nodeAssigne[ nodeID[i] ] = nodeDim;
612               }
613             }
614             if ( token != "Corners" && aHexoticElement )
615               theMesh->SetMeshElementOnShape( aHexoticElement, shapeID );
616           }
617         }
618         delete [] node;
619         delete [] nodeID;
620         break;
621       }
622       case 8: { // "End"
623         EndOfFile = 1;
624         MESSAGE("End of " << theFile << " file");
625         break;
626       }
627       default: {
628         MESSAGE("Unknown Token: " << token);
629       }
630     }
631   }
632   cout << std::endl;
633
634   // remove nodes in holes
635   if ( nbDomains > 1 )
636   {
637     SMESHDS_SubMesh* subMesh;
638     for ( int i = 1; i <= nbVertices; ++i )
639       if ( HexoticNode[i]->NbInverseElements() == 0 )
640       {
641         subMesh =  HexoticNode[i]->getshapeId() > 0 ? theMesh->MeshElements(HexoticNode[i]->getshapeId() ) : 0;
642         theMesh->RemoveFreeNode( HexoticNode[i], subMesh, /*fromGroups=*/false );
643       }
644   }
645   delete [] tabID;
646   delete [] tabRef;
647   delete [] tabDummy;
648   delete [] tabCorner;
649   delete [] tabEdge;
650   delete [] nodeAssigne;
651   delete [] HexoticNode;
652   return true;
653 }
654
655
656 //=======================================================================
657 //function : readResult
658 //purpose  : Read GMF file in case of a mesh w/o geometry
659 //=======================================================================
660
661 static bool readResult(std::string theFile,
662                        HexoticPlugin_Hexotic*  theAlgo,
663                        SMESH_MesherHelper* theHelper)
664 {
665   SMESHDS_Mesh* theMesh = theHelper->GetMeshDS();
666
667   // ---------------------------------
668   // Read generated elements and nodes
669   // ---------------------------------
670
671   std::string token;
672   const int nbField = 9;
673   int nField, EndOfFile = 0, nbElem = 0, nbRef = 0;
674   int aHexoticNodeID = 0, shapeID;
675   int tabRef[nbField], *nodeAssigne;
676   bool tabDummy[nbField], hasDummy = false;
677   std::map <std::string,int> mapField;
678   SMDS_MeshNode** HexoticNode;
679
680   mapField["MeshVersionFormatted"] = 0; tabRef[0] = 0; tabDummy[0] = false;
681   mapField["Dimension"]            = 1; tabRef[1] = 0; tabDummy[1] = false;
682   mapField["Vertices"]             = 2; tabRef[2] = 3; tabDummy[2] = true;
683   mapField["Corners"]              = 3; tabRef[3] = 1; tabDummy[3] = false;
684   mapField["Edges"]                = 4; tabRef[4] = 2; tabDummy[4] = true;
685   mapField["Ridges"]               = 5; tabRef[5] = 1; tabDummy[5] = false;
686   mapField["Quadrilaterals"]       = 6; tabRef[6] = 4; tabDummy[6] = true;
687   mapField["Hexahedra"]            = 7; tabRef[7] = 8; tabDummy[7] = true;
688   mapField["End"]                  = 8; tabRef[8] = 0; tabDummy[8] = false;
689
690   {
691     // theMesh->Clear(); -- this does not remove imported mesh
692     SMDS_ElemIteratorPtr eIt = theMesh->elementsIterator();
693     while( eIt->more() )
694       theMesh->RemoveFreeElement( eIt->next(), /*sm=*/0 );
695     SMDS_NodeIteratorPtr nIt = theMesh->nodesIterator();
696     while ( nIt->more() )
697       theMesh->RemoveFreeNode( nIt->next(), /*sm=*/0 );
698   }
699
700   int nbVertices = getNbShape(theFile, "Vertices");
701   HexoticNode = new SMDS_MeshNode*[ nbVertices + 1 ];
702   nodeAssigne = new int[ nbVertices + 1 ];
703
704   MESSAGE("Read " << theFile << " file");
705   std::ifstream fileRes(theFile.c_str());
706   ASSERT(fileRes);
707
708   while ( !EndOfFile  )
709   {
710     int dummy;
711     fileRes >> token;
712
713     if (mapField.count(token)) {
714       nField   = mapField[token];
715       nbRef    = tabRef[nField];
716       hasDummy = tabDummy[nField];
717     }
718     else {
719       nField = -1;
720       nbRef = 0;
721     }
722
723     nbElem = 0;
724     if ( nField < (mapField.size() - 1) && nField >= 0 )
725       fileRes >> nbElem;
726
727     switch (nField) {
728     case 0: { // "MeshVersionFormatted"
729       MESSAGE(token << " " << nbElem);
730       break;
731     }
732     case 1: { // "Dimension"
733       MESSAGE("Mesh dimension " << nbElem << "D");
734       break;
735     }
736     case 2: { // "Vertices"
737       MESSAGE("Read " << nbElem << " " << token);
738       int aHexoticID;
739       double coord[3];
740       SMDS_MeshNode * aHexoticNode;
741
742       for ( int iElem = 0; iElem < nbElem; iElem++ ) {
743         if(theAlgo->computeCanceled())
744           {
745             return false;
746           }
747         aHexoticID = iElem + 1;
748         for ( int iCoord = 0; iCoord < 3; iCoord++ )
749           fileRes >> coord[ iCoord ];
750         fileRes >> dummy;
751         aHexoticNode = theMesh->AddNode(coord[0], coord[1], coord[2]);
752         HexoticNode[ aHexoticID ] = aHexoticNode;
753         nodeAssigne[ aHexoticID ] = 0;
754       }
755       break;
756     }
757     case 3: // "Corners"
758     case 4: // "Edges"
759     case 5: // "Ridges"
760     case 6: // "Quadrilaterals"
761     case 7: { // "Hexahedra"
762       MESSAGE("Read " << nbElem << " " << token);
763       std::vector< SMDS_MeshNode* > node( nbRef );
764       std::vector< int >          nodeID( nbRef );
765
766       for ( int iElem = 0; iElem < nbElem; iElem++ )
767       {
768         if(theAlgo->computeCanceled())
769           {
770             return false;
771           }
772         for ( int iRef = 0; iRef < nbRef; iRef++ )
773         {
774           fileRes >> aHexoticNodeID;                          // read nbRef aHexoticNodeID
775           node  [ iRef ] = HexoticNode[ aHexoticNodeID ];
776           nodeID[ iRef ] = aHexoticNodeID;
777         }
778         if ( hasDummy )
779           fileRes >> dummy;
780         switch (nField)
781         {
782         case 4: // "Edges"
783           theHelper->AddEdge( node[0], node[1] ); break;
784         case 6:  // "Quadrilaterals"
785           theMesh->AddFace( node[0], node[1], node[2], node[3] ); break;
786         case 7: // "Hexahedra"
787           theHelper->AddVolume( node[0], node[3], node[2], node[1],
788                                 node[4], node[7], node[6], node[5] ); break;
789         default: continue;
790         }
791         if ( nField == 6 )
792           for ( int iRef = 0; iRef < nbRef; iRef++ )
793             nodeAssigne[ nodeID[ iRef ]] = 1;
794       }
795       break;
796     }
797     case 8: { // "End"
798       EndOfFile = 1;
799       MESSAGE("End of " << theFile << " file");
800       break;
801     }
802     default: {
803       MESSAGE("Unknown Token: " << token);
804     }
805     }
806   }
807   cout << std::endl;
808
809   shapeID = theHelper->GetSubShapeID();
810   for ( int i = 0; i < nbVertices; ++i )
811     if ( !nodeAssigne[ i+1 ])
812       theMesh->SetNodeInVolume( HexoticNode[ i+1 ], shapeID );
813
814   delete [] HexoticNode;
815   delete [] nodeAssigne;
816   return true;
817 }
818
819 //=============================================================================
820 /*!
821  * Pass parameters to Hexotic
822  */
823 //=============================================================================
824
825 void HexoticPlugin_Hexotic::SetParameters(const HexoticPlugin_Hypothesis* hyp) {
826
827   MESSAGE("HexoticPlugin_Hexotic::SetParameters");
828   if (hyp) {
829     _hexesMinLevel = hyp->GetHexesMinLevel();
830     _hexesMaxLevel = hyp->GetHexesMaxLevel();
831     _hexesMinSize = hyp->GetMinSize();
832     _hexesMaxSize = hyp->GetMaxSize();
833     _hexoticIgnoreRidges = hyp->GetHexoticIgnoreRidges();
834     _hexoticInvalidElements = hyp->GetHexoticInvalidElements();
835     _hexoticSharpAngleThreshold = hyp->GetHexoticSharpAngleThreshold();
836     _hexoticNbProc = hyp->GetHexoticNbProc();
837     _hexoticWorkingDirectory = hyp->GetHexoticWorkingDirectory();
838     _hexoticVerbosity = hyp->GetHexoticVerbosity();
839     _hexoticMaxMemory = hyp->GetHexoticMaxMemory();
840     _hexoticSdMode = hyp->GetHexoticSdMode();
841     _sizeMaps = hyp->GetSizeMaps();
842   }
843   else {
844     cout << std::endl;
845     cout << "WARNING : The Hexotic default parameters are taken into account" << std::endl;
846     cout << "=======" << std::endl;
847     _hexesMinLevel = hyp->GetDefaultHexesMinLevel();
848     _hexesMaxLevel = hyp->GetDefaultHexesMaxLevel();
849     _hexesMinSize = hyp->GetDefaultMinSize();
850     _hexesMaxSize = hyp->GetDefaultMaxSize();
851     _hexoticIgnoreRidges = hyp->GetDefaultHexoticIgnoreRidges();
852     _hexoticInvalidElements = hyp->GetDefaultHexoticInvalidElements();
853     _hexoticSharpAngleThreshold = hyp->GetDefaultHexoticSharpAngleThreshold();
854     _hexoticNbProc = hyp->GetDefaultHexoticNbProc();
855     _hexoticWorkingDirectory = hyp->GetDefaultHexoticWorkingDirectory();
856     _hexoticVerbosity = hyp->GetDefaultHexoticVerbosity();
857     _hexoticMaxMemory = hyp->GetDefaultHexoticMaxMemory();
858     _hexoticSdMode = hyp->GetDefaultHexoticSdMode();
859     _sizeMaps = hyp->GetDefaultHexoticSizeMaps();
860   }
861 }
862
863 //=======================================================================
864 //function : getTmpDir
865 //purpose  :
866 //=======================================================================
867
868 static TCollection_AsciiString getTmpDir()
869 {
870   TCollection_AsciiString aTmpDir;
871
872   char *Tmp_dir = getenv("SALOME_TMP_DIR");
873 #ifdef WIN32
874   if(Tmp_dir == NULL) {
875     Tmp_dir = getenv("TEMP");
876     if( Tmp_dir== NULL )
877       Tmp_dir = getenv("TMP");
878   }
879 #endif
880
881   if(Tmp_dir != NULL) {
882     aTmpDir = Tmp_dir;
883 #ifdef WIN32
884     if(aTmpDir.Value(aTmpDir.Length()) != '\\') aTmpDir+='\\';
885 #else
886     if(aTmpDir.Value(aTmpDir.Length()) != '/') aTmpDir+='/';
887 #endif
888   }
889   else {
890 #ifdef WIN32
891     aTmpDir = TCollection_AsciiString("C:\\");
892 #else
893     aTmpDir = TCollection_AsciiString("/tmp/");
894 #endif
895   }
896   return aTmpDir;
897 }
898
899 //=======================================================================
900 //function : getSuffix
901 //purpose  : Returns a suffix that will be unique for the current process
902 //=======================================================================
903
904 static TCollection_AsciiString getSuffix()
905 {
906   TCollection_AsciiString aSuffix = "";
907   aSuffix += "_";
908   aSuffix += getenv("USER");
909   aSuffix += "_";
910   aSuffix += Kernel_Utils::GetHostname().c_str();
911   aSuffix += "_";
912 #ifndef WIN32
913   aSuffix += getpid();
914 #else
915   aSuffix += _getpid();
916 #endif
917
918   return aSuffix;
919 }
920
921 //================================================================================
922 /*!
923  * \brief Returns a command to run Hexotic mesher
924  */
925 //================================================================================
926
927 std::string HexoticPlugin_Hexotic::getHexoticCommand(const TCollection_AsciiString& Hexotic_In,
928                                                      const TCollection_AsciiString& Hexotic_Out,
929                                                      const TCollection_AsciiString& Hexotic_SizeMap_Prefix) const
930 {
931   cout << std::endl;
932   cout << "Hexotic execution..." << std::endl;
933   cout << _name << " parameters :" << std::endl;
934   cout << "    " << _name << " Verbosity = " << _hexoticVerbosity << std::endl;
935   cout << "    " << _name << " Max Memory = " << _hexoticMaxMemory << std::endl;
936   cout << "    " << _name << " Segments Min Level = " << _hexesMinLevel << std::endl;
937   cout << "    " << _name << " Segments Max Level = " << _hexesMaxLevel << std::endl;
938   cout << "    " << _name << " Segments Min Size = " << _hexesMinSize << std::endl;
939   cout << "    " << _name << " Segments Max Size = " << _hexesMaxSize << std::endl;
940   cout << "    " << "Hexotic can ignore ridges : " << (_hexoticIgnoreRidges ? "yes":"no") << std::endl;
941   cout << "    " << "Hexotic authorize invalide elements : " << ( _hexoticInvalidElements ? "yes":"no") << std::endl;
942   cout << "    " << _name << " Sharp angle threshold = " << _hexoticSharpAngleThreshold << " degrees" << std::endl;
943   cout << "    " << _name << " Number of threads = " << _hexoticNbProc << std::endl;
944   cout << "    " << _name << " Working directory = \"" << _hexoticWorkingDirectory << "\"" << std::endl;
945   cout << "    " << _name << " Sub. Dom mode = " << _hexoticSdMode << std::endl;
946
947   TCollection_AsciiString run_Hexotic( "mg-hexa.exe" );
948
949   TCollection_AsciiString minl = " --min_level ", maxl = " --max_level ", angle = " --ridge_angle ";
950   TCollection_AsciiString mins = " --min_size ", maxs = " --max_size ";
951   TCollection_AsciiString in   = " --in ",   out  = " --out ";
952   TCollection_AsciiString sizeMap = " --read_sizemap ";
953   TCollection_AsciiString ignoreRidges = " --compute_ridges no ", invalideElements = " --allow_invalid_elements yes ";
954   TCollection_AsciiString subdom = " --components ";
955   TCollection_AsciiString proc = " --max_number_of_threads ";
956   TCollection_AsciiString verb = " --verbose ";
957   TCollection_AsciiString maxmem = " --max_memory ";
958
959   TCollection_AsciiString minLevel, maxLevel, minSize, maxSize, sharpAngle, mode, nbproc, verbosity, maxMemory;
960   minLevel = _hexesMinLevel;
961   maxLevel = _hexesMaxLevel;
962   minSize = _hexesMinSize;
963   maxSize = _hexesMaxSize;
964   sharpAngle = _hexoticSharpAngleThreshold;
965   // Mode translation for mg-tetra 1.1
966   switch ( _hexoticSdMode )
967   {
968     case 1:
969       mode = "outside_skin_only";
970       break;
971     case 2:
972       mode = "outside_components";
973       break;
974     case 3:
975       mode = "all";
976       break;
977     case 4:
978       mode = "all --manifold_geometry no";
979       break;
980   }
981   nbproc = _hexoticNbProc;
982   verbosity = _hexoticVerbosity;
983   maxMemory = _hexoticMaxMemory;
984
985   if (_hexoticIgnoreRidges)
986     run_Hexotic +=  ignoreRidges;
987
988   if (_hexoticInvalidElements)
989     run_Hexotic +=  invalideElements;
990
991   if (_hexesMinSize > 0)
992     run_Hexotic +=  mins + minSize;
993
994   if (_hexesMaxSize > 0)
995     run_Hexotic +=  maxs + maxSize;
996
997   if (_hexesMinLevel > 0)
998     run_Hexotic +=  minl + minLevel;
999
1000   if (_hexesMaxLevel > 0)
1001     run_Hexotic +=  maxl + maxLevel;
1002
1003   if (_hexoticSharpAngleThreshold > 0)
1004     run_Hexotic +=  angle + sharpAngle;
1005   
1006   if (_sizeMaps.begin() != _sizeMaps.end())
1007     run_Hexotic += sizeMap + Hexotic_SizeMap_Prefix;
1008
1009   run_Hexotic += in + Hexotic_In + out + Hexotic_Out;
1010   run_Hexotic += subdom + mode;
1011   run_Hexotic += proc + nbproc;
1012   run_Hexotic += verb + verbosity;
1013   run_Hexotic += maxmem + maxMemory;
1014
1015   return run_Hexotic.ToCString();
1016 }
1017
1018 // TODO : this is a duplication of some code found in BLSURFPlugin_BLSURF find a proper
1019 // way to share it
1020 TopoDS_Shape HexoticPlugin_Hexotic::entryToShape(std::string entry)
1021 {
1022   MESSAGE("HexoticPlugin_Hexotic::entryToShape "<<entry );
1023   GEOM::GEOM_Object_var aGeomObj;
1024   TopoDS_Shape S = TopoDS_Shape();
1025   SALOMEDS::SObject_var aSObj = myStudy->FindObjectID( entry.c_str() );
1026   if (!aSObj->_is_nil()) {
1027     CORBA::Object_var obj = aSObj->GetObject();
1028     aGeomObj = GEOM::GEOM_Object::_narrow(obj);
1029     aSObj->UnRegister();
1030   }
1031   if ( !aGeomObj->_is_nil() )
1032     S = smeshGen_i->GeomObjectToShape( aGeomObj.in() );
1033   return S;
1034 }
1035
1036 //================================================================================
1037 /*!
1038  * \brief Produces a .mesh file with the size maps informations to give to Hexotic
1039  */
1040 //================================================================================
1041 std::vector<std::string> HexoticPlugin_Hexotic::writeSizeMapFile( std::string sizeMapPrefix )
1042 {
1043   HexoticPlugin_Hypothesis::THexoticSizeMaps::iterator it ;
1044   
1045   // The GMF driver will be used to write the size map file
1046   DriverGMF_Write aWriter;
1047   aWriter.SetSizeMapPrefix( sizeMapPrefix ); 
1048   
1049   std::vector<Control_Pnt> points;
1050   // Iterate on the size maps
1051   for (it=_sizeMaps.begin(); it!=_sizeMaps.end(); it++)
1052   {
1053     // Step 1 : Get the GEOM object entry and the size 
1054     // from the _sizeMaps infos
1055     std::string anEntry = it->first;
1056     double aLocalSize = it->second;
1057     TopoDS_Shape aShape = entryToShape( anEntry );
1058     
1059     // Step 2 : Create the points
1060     createControlPoints( aShape, aLocalSize, points );
1061   }
1062   // Write the .mesh size map file
1063   aWriter.PerformSizeMap(points); 
1064   return aWriter.GetSizeMapFiles();
1065 }
1066
1067 //================================================================================
1068 /*!
1069  * \brief Fills a vector of points from which a size map input file can be written
1070  */
1071 //================================================================================
1072 void HexoticPlugin_Hexotic::createControlPoints( const TopoDS_Shape& aShape, 
1073                                                  const double& theSize, 
1074                                                  std::vector<Control_Pnt>& thePoints )
1075
1076   if ( aShape.ShapeType() == TopAbs_VERTEX )
1077   {
1078     gp_Pnt aPnt = BRep_Tool::Pnt( TopoDS::Vertex(aShape) );
1079     Control_Pnt aControl_Pnt( aPnt, theSize );
1080     thePoints.push_back( aControl_Pnt );
1081   }
1082   if ( aShape.ShapeType() == TopAbs_EDGE )
1083   {
1084     createPointsSampleFromEdge( aShape, theSize, thePoints );  
1085   }
1086   else if ( aShape.ShapeType() == TopAbs_WIRE )
1087   {
1088     TopExp_Explorer Ex;
1089     for (Ex.Init(aShape,TopAbs_EDGE); Ex.More(); Ex.Next()) 
1090     {
1091       createPointsSampleFromEdge( Ex.Current(), theSize, thePoints );
1092     } 
1093   }
1094   else if ( aShape.ShapeType() ==  TopAbs_FACE )
1095   {
1096     createPointsSampleFromFace( aShape, theSize, thePoints ); 
1097   }
1098   else if ( aShape.ShapeType() ==  TopAbs_SOLID )
1099   {
1100     createPointsSampleFromSolid( aShape, theSize, thePoints ); 
1101   }
1102   else if ( aShape.ShapeType() == TopAbs_COMPOUND )
1103   {
1104     TopoDS_Iterator it( aShape );
1105     for(; it.More(); it.Next())
1106     {
1107       createControlPoints( it.Value(), theSize, thePoints );
1108     }
1109   }
1110 }
1111
1112 //================================================================================
1113 /*!
1114  * \brief Fills a vector of points with point samples approximately 
1115  * \brief spaced with a given size
1116  */
1117 //================================================================================
1118 void HexoticPlugin_Hexotic::createPointsSampleFromEdge( const TopoDS_Shape& aShape, 
1119                                                         const double& theSize, 
1120                                                         std::vector<Control_Pnt>& thePoints )
1121 {
1122   double step = theSize;
1123   double first, last;  
1124   Handle( Geom_Curve ) aCurve = BRep_Tool::Curve( TopoDS::Edge( aShape ), first, last );
1125   GeomAdaptor_Curve C ( aCurve );
1126   GCPnts_UniformAbscissa DiscretisationAlgo(C, step , first, last, Precision::Confusion());
1127   int nbPoints = DiscretisationAlgo.NbPoints();
1128   
1129   for ( int i = 1; i <= nbPoints; i++ )
1130   {
1131     double param = DiscretisationAlgo.Parameter( i );
1132     Control_Pnt aPnt;
1133     aCurve->D0( param, aPnt );
1134     aPnt.SetSize(theSize);
1135     thePoints.push_back( aPnt );
1136   }  
1137 }
1138
1139 //================================================================================
1140 /*!
1141  * \brief Fills a vector of points with point samples approximately 
1142  * \brief spaced with a given size
1143  */
1144 //================================================================================
1145 void HexoticPlugin_Hexotic::createPointsSampleFromFace( const TopoDS_Shape& aShape, 
1146                                                         const double& theSize, 
1147                                                         std::vector<Control_Pnt>& thePoints )
1148 {
1149   BRepMesh_IncrementalMesh M(aShape, 0.01, Standard_True);
1150   TopLoc_Location aLocation;
1151   TopoDS_Face aFace = TopoDS::Face(aShape);
1152
1153   // Triangulate the face
1154   Handle(Poly_Triangulation) aTri = BRep_Tool::Triangulation (aFace, aLocation);
1155   
1156   // Get the transformation associated to the face location
1157   gp_Trsf aTrsf = aLocation.Transformation();
1158   
1159   // Get triangles
1160   int nbTriangles = aTri->NbTriangles();
1161   Poly_Array1OfTriangle triangles(1,nbTriangles);
1162   triangles=aTri->Triangles();
1163   
1164   // GetNodes
1165   int nbNodes = aTri->NbNodes();
1166   TColgp_Array1OfPnt nodes(1,nbNodes);
1167   nodes = aTri->Nodes();
1168
1169   // Iterate on triangles and subdivide them
1170   for(int i=1; i<=nbTriangles; i++)
1171   {
1172      Poly_Triangle aTriangle = triangles.Value(i);
1173      gp_Pnt p1 = nodes.Value(aTriangle.Value(1));
1174      gp_Pnt p2 = nodes.Value(aTriangle.Value(2));
1175      gp_Pnt p3 = nodes.Value(aTriangle.Value(3));
1176      
1177      p1.Transform(aTrsf);
1178      p2.Transform(aTrsf);
1179      p3.Transform(aTrsf);
1180      
1181      subdivideTriangle(p1, p2, p3, theSize, thePoints);  
1182   }
1183 }
1184
1185 //================================================================================
1186 /*!
1187  * \brief Fills a vector of points with point samples approximately 
1188  * \brief spaced with a given size
1189  */
1190 //================================================================================
1191 void HexoticPlugin_Hexotic::createPointsSampleFromSolid( const TopoDS_Shape& aShape, 
1192                                                          const double& theSize, 
1193                                                          std::vector<Control_Pnt>& thePoints )
1194 {
1195   // Compute the bounding box
1196   double Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
1197   Bnd_Box B;               
1198   BRepBndLib::Add(aShape, B);
1199   B.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
1200   
1201   // Create the points
1202   double step = theSize;
1203   
1204   for ( double x=Xmin; x-Xmax<Precision::Confusion(); x=x+step )
1205   {
1206     for ( double y=Ymin; y-Ymax<Precision::Confusion(); y=y+step )
1207     {
1208       // Step1 : generate the Zmin -> Zmax line
1209       gp_Pnt startPnt(x, y, Zmin);
1210       gp_Pnt endPnt(x, y, Zmax);
1211       gp_Vec aVec(startPnt, endPnt);
1212       gp_Lin aLine(startPnt, aVec);
1213       double endParam = Zmax - Zmin;
1214       
1215       // Step2 : for each face of aShape:
1216       std::set<double> intersections;
1217       std::set<double>::iterator it = intersections.begin();
1218       
1219       TopExp_Explorer Ex;
1220       for (Ex.Init(aShape,TopAbs_FACE); Ex.More(); Ex.Next()) 
1221       { 
1222         // check if there is an intersection
1223         IntCurvesFace_Intersector anIntersector(TopoDS::Face(Ex.Current()), Precision::Confusion());
1224         anIntersector.Perform(aLine, 0, endParam);
1225         
1226         // get the intersection's parameter and store it
1227         int nbPoints = anIntersector.NbPnt();
1228         for(int i = 0 ; i < nbPoints ; i++ )
1229         {
1230           it = intersections.insert( it, anIntersector.WParameter(i+1) );
1231         }
1232       }
1233       // Step3 : go through the line chunk by chunk 
1234       if ( intersections.begin() != intersections.end() )
1235       {
1236         std::set<double>::iterator intersectionsIterator=intersections.begin();
1237         double first = *intersectionsIterator;
1238         intersectionsIterator++;
1239         bool innerPoints = true; 
1240         for ( ; intersectionsIterator!=intersections.end() ; intersectionsIterator++ )
1241         {
1242           double second = *intersectionsIterator;
1243           if ( innerPoints )
1244           {
1245             // If the last chunk was outside of the shape or this is the first chunk
1246             // add the points in the range [first, second] to the points vector
1247             double localStep = (second -first) / ceil( (second - first) / step );
1248             for ( double z = Zmin + first; z < Zmin + second; z = z + localStep )
1249             {
1250               thePoints.push_back(Control_Pnt( x, y, z, theSize ));
1251             }
1252             thePoints.push_back(Control_Pnt( x, y, Zmin + second, theSize ));
1253           }
1254           first = second;
1255           innerPoints = !innerPoints;
1256         }
1257       }
1258     }
1259   }
1260 }
1261
1262 //================================================================================
1263 /*!
1264  * \brief Subdivides a triangle until it reaches a certain size (recursive function)
1265  */
1266 //================================================================================
1267 void HexoticPlugin_Hexotic::subdivideTriangle( const gp_Pnt& p1, 
1268                                                const gp_Pnt& p2, 
1269                                                const gp_Pnt& p3, 
1270                                                const double& theSize, 
1271                                                std::vector<Control_Pnt>& thePoints)
1272 {
1273   // Size threshold to stop subdividing
1274   // This value ensures that two control points are distant no more than 2*theSize
1275   // as shown below
1276   //
1277   // The greater distance D of the mass center M to each Edge is 1/3 * Median 
1278   // and Median < sqrt(3/4) * a  where a is the greater side (by using Apollonius' thorem). 
1279   // So D < 1/3 * sqrt(3/4) * a and if a < sqrt(3) * S then D < S/2
1280   // and the distance between two mass centers of two neighbouring triangles 
1281   // sharing an edge is < 2 * 1/2 * S = S
1282   // If the traingles share a Vertex and no Edge the distance of the mass centers 
1283   // to the Vertices is 2*D < S so the mass centers are distant of less than 2*S 
1284   
1285   double threshold = sqrt( 3. ) * theSize;
1286   
1287   if ( (p1.Distance(p2) > threshold ||
1288         p2.Distance(p3) > threshold ||
1289         p3.Distance(p1) > threshold))
1290   { 
1291     std::vector<gp_Pnt> midPoints = computePointsForSplitting(p1, p2, p3);
1292  
1293     subdivideTriangle( midPoints[0], midPoints[1], midPoints[2], theSize, thePoints );
1294     subdivideTriangle( midPoints[0], p2, midPoints[1], theSize, thePoints );
1295     subdivideTriangle( midPoints[2], midPoints[1], p3, theSize, thePoints );
1296     subdivideTriangle( p1, midPoints[0], midPoints[2], theSize, thePoints );
1297   }
1298   else
1299   {
1300     double x = (p1.X() + p2.X() + p3.X()) / 3 ;
1301     double y = (p1.Y() + p2.Y() + p3.Y()) / 3 ;
1302     double z = (p1.Z() + p2.Z() + p3.Z()) / 3 ;
1303     
1304     Control_Pnt massCenter( x ,y ,z, theSize );
1305     thePoints.push_back( massCenter );
1306   }
1307 }
1308
1309 //================================================================================
1310 /*!
1311  * \brief Returns the appropriate points for splitting a triangle
1312  * \brief the tangency points of the incircle are used in order to have mostly
1313  * \brief well-shaped sub-triangles
1314  */
1315 //================================================================================
1316 std::vector<gp_Pnt> HexoticPlugin_Hexotic::computePointsForSplitting( const gp_Pnt& p1, 
1317                                                                       const gp_Pnt& p2, 
1318                                                                       const gp_Pnt& p3 )
1319 {
1320   std::vector<gp_Pnt> midPoints;
1321   //Change coordinates
1322   gp_Trsf Trsf_1;            // Identity transformation
1323   gp_Ax3 reference_system(gp::Origin(), gp::DZ(), gp::DX());   // OXY
1324  
1325   gp_Vec Vx(p1, p3);
1326   gp_Vec Vaux(p1, p2);
1327   gp_Dir Dx(Vx);
1328   gp_Dir Daux(Vaux);
1329   gp_Dir Dz = Dx.Crossed(Daux);
1330   gp_Ax3 current_system(p1, Dz, Dx);
1331   
1332   Trsf_1.SetTransformation( reference_system, current_system );
1333   
1334   gp_Pnt A = p1.Transformed(Trsf_1);
1335   gp_Pnt B = p2.Transformed(Trsf_1);
1336   gp_Pnt C = p3.Transformed(Trsf_1);
1337   
1338   double a =  B.Distance(C) ;
1339   double b =  A.Distance(C) ;
1340   double c =  B.Distance(A) ;
1341   
1342   // Incenter coordinates
1343   // see http://mathworld.wolfram.com/Incenter.html
1344   double Xi = ( b*B.X() + c*C.X() ) / ( a + b + c );
1345   double Yi = ( b*B.Y() ) / ( a + b + c );
1346   gp_Pnt Center(Xi, Yi, 0);
1347   
1348   // Calculate the tangency points of the incircle
1349   gp_Pnt T1 = tangencyPoint( A, B, Center);
1350   gp_Pnt T2 = tangencyPoint( B, C, Center);
1351   gp_Pnt T3 = tangencyPoint( C, A, Center);
1352   
1353   gp_Pnt p1_2 = T1.Transformed(Trsf_1.Inverted());
1354   gp_Pnt p2_3 = T2.Transformed(Trsf_1.Inverted());
1355   gp_Pnt p3_1 = T3.Transformed(Trsf_1.Inverted());
1356
1357   midPoints.push_back(p1_2);
1358   midPoints.push_back(p2_3);
1359   midPoints.push_back(p3_1);
1360   
1361   return midPoints;
1362 }
1363
1364 //================================================================================
1365 /*!
1366  * \brief Computes the tangency points of the circle of center Center with
1367  * \brief the straight line (p1 p2)
1368  */
1369 //================================================================================
1370 gp_Pnt HexoticPlugin_Hexotic::tangencyPoint(const gp_Pnt& p1,
1371                                             const gp_Pnt& p2,
1372                                             const gp_Pnt& Center)
1373 {
1374   double Xt = 0;
1375   double Yt = 0;
1376   
1377   // The tangency point is the intersection of the straight line (p1 p2)
1378   // and the straight line (Center T) which is orthogonal to (p1 p2)
1379   if ( fabs(p1.X() - p2.X()) <= Precision::Confusion() )
1380   {
1381     Xt=p1.X();     // T is on (p1 p2)
1382     Yt=Center.Y(); // (Center T) is orthogonal to (p1 p2)
1383   }
1384   else if ( fabs(p1.Y() - p2.Y()) <= Precision::Confusion() )
1385   {
1386     Yt=p1.Y();     // T is on (p1 p2) 
1387     Xt=Center.X(); // (Center T) is orthogonal to (p1 p2)
1388   }
1389   else
1390   {
1391     // First straight line coefficients (equation y=a*x+b)
1392     double a = (p2.Y() - p1.Y()) / (p2.X() - p1.X())  ;
1393     double b = p1.Y() - a*p1.X();         // p1 is on this straight line
1394     
1395     // Second straight line coefficients (equation y=c*x+d)
1396     double c = -1 / a;                    // The 2 lines are orthogonal
1397     double d = Center.Y() - c*Center.X(); // Center is on this straight line
1398     
1399     Xt = (d - b) / (a - c);
1400     Yt = a*Xt + b;
1401   }
1402   
1403   return gp_Pnt( Xt, Yt, 0 );
1404 }
1405
1406 //=============================================================================
1407 /*!
1408  * Here we are going to use the Hexotic mesher
1409  */
1410 //=============================================================================
1411
1412 bool HexoticPlugin_Hexotic::Compute(SMESH_Mesh&          aMesh,
1413                                      const TopoDS_Shape& aShape)
1414 {
1415   _compute_canceled = false;
1416   bool Ok = true;
1417   SMESHDS_Mesh* meshDS = aMesh.GetMeshDS();
1418   TCollection_AsciiString hexahedraMessage;
1419
1420   if (_iShape == 0 && _nbShape == 0) {
1421     _nbShape = countShape( meshDS, TopAbs_SOLID );  // we count the number of shapes
1422   }
1423
1424   // to prevent from displaying error message after computing,
1425   // SetIsAlwaysComputed( true ) to empty sub-meshes
1426   vector< SMESH_subMesh* > subMeshesAlwaysComp;
1427   for ( int i = 0; i < _nbShape; ++i )
1428     if ( SMESH_subMesh* sm = aMesh.GetSubMeshContaining( aShape ))
1429     {
1430       SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(/*includeSelf=*/true,
1431                                                                /*complexShapeFirst=*/false);
1432       while ( smIt->more() )
1433       {
1434         sm = smIt->next();
1435         if ( !sm->IsMeshComputed() )
1436         {
1437           sm->SetIsAlwaysComputed( true );
1438           subMeshesAlwaysComp.push_back( sm );
1439         }
1440       }
1441     }
1442
1443   _iShape++;
1444
1445   if (_iShape == _nbShape ) {
1446
1447     // create bounding box for each shape of the compound
1448
1449     int iShape = 0;
1450     TopoDS_Shape *tabShape;
1451     double **tabBox;
1452
1453     tabShape = new TopoDS_Shape[_nbShape];
1454     tabBox   = new double*[_nbShape];
1455     for (int i=0; i<_nbShape; i++)
1456       tabBox[i] = new double[6];
1457     double Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
1458
1459     TopExp_Explorer expBox (meshDS->ShapeToMesh(), TopAbs_SOLID);
1460     for (; expBox.More(); expBox.Next()) {
1461       tabShape[iShape] = expBox.Current();
1462       Bnd_Box BoundingBox;
1463       BRepBndLib::Add(expBox.Current(), BoundingBox);
1464       BoundingBox.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
1465       tabBox[iShape][0] = Xmin; tabBox[iShape][1] = Xmax;
1466       tabBox[iShape][2] = Ymin; tabBox[iShape][3] = Ymax;
1467       tabBox[iShape][4] = Zmin; tabBox[iShape][5] = Zmax;
1468       iShape++;
1469     }
1470
1471     SetParameters(_hypothesis);
1472
1473 //     TCollection_AsciiString aTmpDir = getTmpDir();
1474     TCollection_AsciiString aTmpDir = TCollection_AsciiString(_hexoticWorkingDirectory.c_str());
1475 #ifdef WIN32
1476     if ( aTmpDir.Value(aTmpDir.Length()) != '\\' ) aTmpDir += '\\';
1477 #else
1478     if ( aTmpDir.Value(aTmpDir.Length()) != '/' ) aTmpDir += '/';
1479 #endif
1480     TCollection_AsciiString Hexotic_In(""), Hexotic_Out, Hexotic_SizeMap_Prefix;
1481     TCollection_AsciiString modeFile_In( "chmod 666 " ), modeFile_Out( "chmod 666 " );
1482     TCollection_AsciiString aLogFileName = aTmpDir + "Hexotic"+getSuffix()+".log";    // log
1483
1484     std::map <int,int> aSmdsToHexoticIdMap;
1485     std::map <int,const SMDS_MeshNode*> aHexoticIdToNodeMap;
1486
1487     Hexotic_Out = aTmpDir + "Hexotic"+getSuffix()+"_Out.mesh";
1488 #ifdef WITH_BLSURFPLUGIN
1489     bool defaultInputFile = true;
1490     if (_blsurfHypo && !_blsurfHypo->GetQuadAllowed()) {
1491       Hexotic_In = TCollection_AsciiString(_blsurfHypo->GetGMFFile().c_str());
1492       if (Hexotic_In != "")
1493         defaultInputFile = false;
1494     }
1495     if (defaultInputFile) {
1496 #endif
1497       Hexotic_In  = aTmpDir + "Hexotic"+getSuffix()+"_In.mesh";
1498       removeHexoticFiles(Hexotic_In, Hexotic_Out);
1499       cout << std::endl;
1500       cout << "Creating Hexotic input mesh file : " << Hexotic_In << std::endl;
1501       aMesh.ExportGMF(Hexotic_In.ToCString(), meshDS, true);
1502 #ifdef WITH_BLSURFPLUGIN
1503     }
1504     else {
1505       removeFile( Hexotic_Out );
1506     }
1507 #endif
1508     
1509     Hexotic_SizeMap_Prefix = aTmpDir + "Hexotic_SizeMap" + getSuffix();
1510     std::vector<std::string> sizeMapFiles = writeSizeMapFile( Hexotic_SizeMap_Prefix.ToCString() );
1511     
1512     std::string run_Hexotic = getHexoticCommand(Hexotic_In, Hexotic_Out, Hexotic_SizeMap_Prefix);
1513     run_Hexotic += std::string(" 1> ") + aLogFileName.ToCString();  // dump into file
1514
1515     cout << std::endl;
1516     cout << "Hexotic command : " << run_Hexotic << std::endl;
1517
1518 #ifndef WIN32    
1519     modeFile_In += Hexotic_In;
1520     system( modeFile_In.ToCString() );
1521 #endif
1522     aSmdsToHexoticIdMap.clear();
1523     aHexoticIdToNodeMap.clear();
1524
1525     MESSAGE("HexoticPlugin_Hexotic::Compute");
1526
1527     int status = system( run_Hexotic.data() );
1528
1529     // --------------
1530     // read a result
1531     // --------------
1532
1533     std::ifstream fileRes( Hexotic_Out.ToCString() );
1534     modeFile_Out += Hexotic_Out;
1535     system( modeFile_Out.ToCString() );
1536     if ( ! fileRes.fail() ) {
1537       Ok = readResult( Hexotic_Out.ToCString(),
1538                        this,
1539                        meshDS, _nbShape, tabShape, tabBox );
1540       if(Ok) {
1541 /*********************
1542 // TODO: Detect and remove elements in holes in case of sd mode = 4
1543       // Remove previous nodes and elements
1544       SMDS_ElemIteratorPtr itElement = meshDS->elementsIterator();
1545       SMDS_NodeIteratorPtr itNode = meshDS->nodesIterator();
1546     
1547       while ( itElement->more() )
1548         meshDS->RemoveElement( itElement->next() );
1549       while ( itNode->more() )
1550         meshDS->RemoveNode( itNode->next() );
1551   
1552       SMESH_ComputeErrorPtr myError = aMesh.GMFToMesh(Hexotic_Out.ToCString());
1553       if (myError)
1554 */
1555         hexahedraMessage = "success";
1556         #ifndef _DEBUG_
1557         removeFile(Hexotic_Out);
1558         removeFile(Hexotic_In);
1559         removeFile(aLogFileName);
1560         for( int i=0; i<sizeMapFiles.size(); i++)
1561         {
1562           removeFile( TCollection_AsciiString( sizeMapFiles[i].c_str() ) );
1563         }
1564         #endif
1565       }
1566       else {
1567         hexahedraMessage = "failed"; 
1568       }
1569     }
1570     else {
1571       hexahedraMessage = "failed";
1572       cout << "Problem with Hexotic output file " << Hexotic_Out.ToCString() << std::endl;
1573       Ok = false;
1574       // analyse log file
1575       SMESH_File logFile( aLogFileName.ToCString() );
1576       if ( !logFile.eof() )
1577       {
1578         char msgLic[] = " Dlim ";
1579         const char* fileBeg = logFile.getPos(), *fileEnd = fileBeg + logFile.size();
1580         if ( std::search( fileBeg, fileEnd, msgLic, msgLic+strlen(msgLic)) != fileEnd )
1581           error("Licence problems.");
1582       }
1583 #ifndef WIN32
1584       if ( status > 0 && WEXITSTATUS(status) == 127 )
1585         error("hexotic: command not found");
1586 #else
1587       int err = errno;
1588       if ( status == 0 && err == ENOENT ) {
1589         error("hexotic: command not found");
1590       }
1591 #endif
1592     }
1593     cout << "Hexahedra meshing " << hexahedraMessage << std::endl;
1594     cout << std::endl;
1595
1596     // restore "always computed" flag of sub-meshes (0022127)
1597     for  ( size_t iSM = 0; iSM < subMeshesAlwaysComp.size(); ++iSM )
1598       subMeshesAlwaysComp[ iSM ]->SetIsAlwaysComputed( false );
1599
1600     delete [] tabShape;
1601     for (int i=0; i<_nbShape; i++)
1602       delete [] tabBox[i];
1603     delete [] tabBox;
1604     _nbShape = 0;
1605     _iShape  = 0;
1606   }
1607   if(_compute_canceled)
1608     return error(SMESH_Comment("interruption initiated by user"));
1609   return Ok;
1610 }
1611
1612 //=============================================================================
1613 /*!
1614  * \brief Computes mesh without geometry
1615  *  \param aMesh - the mesh
1616  *  \param aHelper - helper that must be used for adding elements to \aaMesh
1617  *  \retval bool - is a success
1618  *
1619  * The method is called if ( !aMesh->HasShapeToMesh() )
1620  */
1621 //=============================================================================
1622
1623 bool HexoticPlugin_Hexotic::Compute(SMESH_Mesh & aMesh, SMESH_MesherHelper* aHelper)
1624 {
1625   _compute_canceled = false;
1626 /*
1627   SMESH_ComputeErrorPtr myError = SMESH_ComputeError::New();
1628 */
1629   bool Ok = true;
1630   TCollection_AsciiString hexahedraMessage;
1631
1632   SetParameters(_hypothesis);
1633
1634   TCollection_AsciiString aTmpDir = getTmpDir();
1635   TCollection_AsciiString Hexotic_In, Hexotic_Out, Hexotic_SizeMap_Prefix;
1636   TCollection_AsciiString modeFile_In( "chmod 666 " ), modeFile_Out( "chmod 666 " );
1637   TCollection_AsciiString aLogFileName = aTmpDir + "Hexotic"+getSuffix()+".log";    // log
1638
1639   std::map <int,int> aSmdsToHexoticIdMap;
1640   std::map <int,const SMDS_MeshNode*> aHexoticIdToNodeMap;
1641
1642   Hexotic_In  = aTmpDir + "Hexotic"+getSuffix()+"_In.mesh";
1643   Hexotic_Out = aTmpDir + "Hexotic"+getSuffix()+"_Out.mesh";
1644   Hexotic_SizeMap_Prefix = aTmpDir + "Hexotic_SizeMap";
1645  
1646   
1647   std::vector<std::string> sizeMapFiles = writeSizeMapFile( Hexotic_SizeMap_Prefix.ToCString() );
1648
1649   std::string run_Hexotic = getHexoticCommand(Hexotic_In, Hexotic_Out, Hexotic_SizeMap_Prefix);
1650   run_Hexotic += std::string(" 1> ") + aLogFileName.ToCString();  // dump into file
1651
1652   removeHexoticFiles(Hexotic_In, Hexotic_Out);
1653
1654   cout << std::endl;
1655   cout << "Creating Hexotic input mesh file : " << Hexotic_In << std::endl;
1656   aMesh.ExportGMF(Hexotic_In.ToCString(), aHelper->GetMeshDS());
1657 #ifndef WIN32    
1658   modeFile_In += Hexotic_In;
1659   system( modeFile_In.ToCString() );
1660 #endif
1661   aSmdsToHexoticIdMap.clear();
1662   aHexoticIdToNodeMap.clear();
1663
1664   MESSAGE("HexoticPlugin_Hexotic::Compute");
1665
1666   cout << std::endl;
1667   cout << "Hexotic command : " << run_Hexotic << std::endl;
1668   system( run_Hexotic.data() );
1669
1670   // --------------
1671   // read a result
1672   // --------------
1673
1674   std::ifstream fileRes( Hexotic_Out.ToCString() );
1675   modeFile_Out += Hexotic_Out;
1676   system( modeFile_Out.ToCString() );
1677   if ( ! fileRes.fail() ) {
1678     Ok = readResult( Hexotic_Out.ToCString(),
1679                      this,
1680                      aHelper );
1681     if(Ok)
1682 /*
1683     // Remove previous nodes and elements
1684     SMDS_ElemIteratorPtr itElement = aHelper->GetMeshDS()->elementsIterator();
1685     SMDS_NodeIteratorPtr itNode = aHelper->GetMeshDS()->nodesIterator();
1686     
1687     while ( itElement->more() )
1688       aHelper->GetMeshDS()->RemoveElement( itElement->next() );
1689     while ( itNode->more() )
1690       aHelper->GetMeshDS()->RemoveNode( itNode->next() );
1691
1692     // Import GMF mesh
1693     myError = aMesh.GMFToMesh(Hexotic_Out.ToCString());
1694     
1695     itElement = aHelper->GetMeshDS()->elementsIterator();
1696     itNode = aHelper->GetMeshDS()->nodesIterator();
1697
1698     // Assign nodes and elements to the pseudo shape
1699     while ( itNode->more() )
1700       aHelper->GetMeshDS()->SetNodeInVolume(itNode->next(), 1);
1701     while ( itElement->more() )
1702       aHelper->GetMeshDS()->SetMeshElementOnShape(itElement->next(), 1);
1703
1704     if(myError->IsOK())
1705 */
1706       hexahedraMessage = "success";
1707     else
1708       hexahedraMessage = "failed";
1709   }
1710   else {
1711 /*
1712     myError->myName = COMPERR_EXCEPTION;
1713 */
1714     hexahedraMessage = "failed";
1715     cout << "Problem with Hexotic output file " << Hexotic_Out << std::endl;
1716     // analyse log file
1717     SMESH_File logFile( aLogFileName.ToCString() );
1718     if ( !logFile.eof() )
1719     {
1720       char msgLic[] = " Dlim ";
1721       const char* fileBeg = logFile.getPos(), *fileEnd = fileBeg + logFile.size();
1722       if ( std::search( fileBeg, fileEnd, msgLic, msgLic+strlen(msgLic)) != fileEnd )
1723         return error("Licence problems.");
1724     }
1725     return error(SMESH_Comment("Problem with Hexotic output file ")<<Hexotic_Out);
1726   }
1727   cout << "Hexahedra meshing " << hexahedraMessage << std::endl;
1728   cout << std::endl;
1729
1730   if(_compute_canceled)
1731     return error(SMESH_Comment("interruption initiated by user"));
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 void HexoticPlugin_Hexotic::CancelCompute()
1767 {
1768   _compute_canceled = true;
1769 #ifdef WIN32
1770 #else
1771   TCollection_AsciiString aTmpDir = getTmpDir();
1772   TCollection_AsciiString Hexotic_In = aTmpDir + "Hexotic_In.mesh";
1773   TCollection_AsciiString cmd = TCollection_AsciiString("ps ux | grep ") + Hexotic_In;
1774   cmd += TCollection_AsciiString(" | grep -v grep | awk '{print $2}' | xargs kill -9 > /dev/null 2>&1");
1775   system( cmd.ToCString() );
1776 #endif
1777 }