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