Salome HOME
Merge from V6_main_20120808 08Aug12
[plugins/hexoticplugin.git] / src / HexoticPlugin / HexoticPlugin_Hexotic.cxx
1 // Copyright (C) 2007-2012  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 "SMDS_MeshElement.hxx"
29 #include "SMDS_MeshNode.hxx"
30
31 #include <TopoDS.hxx>
32 #include <TopExp_Explorer.hxx>
33 #include <OSD_File.hxx>
34
35 #include "utilities.h"
36
37 #ifndef WIN32
38 #include <sys/sysinfo.h>
39 #endif
40
41 #ifdef _DEBUG_
42 #define DUMP(txt) \
43 //  cout << txt
44 #else
45 #define DUMP(txt)
46 #endif
47
48 #include <SMESH_Gen.hxx>
49 #include <SMESHDS_Mesh.hxx>
50 #include <SMESHDS_GroupBase.hxx>
51 #include <SMESH_ControlsDef.hxx>
52 #include <SMESH_MesherHelper.hxx>
53 #include "SMESH_HypoFilter.hxx"
54 #include <SMESH_File.hxx>
55
56 #include <list>
57 #include <cstdlib>
58 #include <iostream>
59
60 #include <Standard_ProgramError.hxx>
61
62 #include <BRepClass3d_SolidClassifier.hxx>
63 #include <Bnd_Box.hxx>
64 #include <BRepBndLib.hxx>
65 #include <Precision.hxx>
66
67 #include <BRepBuilderAPI_MakeVertex.hxx>
68 #include <BRep_Tool.hxx>
69 #include <TopTools_MapOfShape.hxx>
70 #include <TopTools_Array1OfShape.hxx>
71 #include <BRepExtrema_DistShapeShape.hxx>
72 #include <TColStd_Array1OfReal.hxx>
73 #include <BRepExtrema_DistShapeShape.hxx>
74
75 static void removeFile( const TCollection_AsciiString& fileName )
76 {
77   try {
78     OSD_File( fileName ).Remove();
79   }
80   catch ( Standard_ProgramError ) {
81     MESSAGE("Can't remove file: " << fileName.ToCString() << " ; file does not exist or permission denied");
82   }
83 }
84
85 //=============================================================================
86 /*!
87  *  
88  */
89 //=============================================================================
90
91 HexoticPlugin_Hexotic::HexoticPlugin_Hexotic(int hypId, int studyId, SMESH_Gen* gen)
92   : SMESH_3D_Algo(hypId, studyId, gen)
93 {
94   MESSAGE("HexoticPlugin_Hexotic::HexoticPlugin_Hexotic");
95   _name = "Hexotic_3D";
96   _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);// 1 bit /shape type
97 //   _onlyUnaryInput = false;
98   _requireShape = false;
99   _iShape=0;
100   _nbShape=0;
101   _hexoticFilesKept=false;
102   _compatibleHypothesis.push_back("Hexotic_Parameters");
103 #ifdef WITH_BLSURFPLUGIN
104   _blsurfHypo = NULL;
105 #endif
106 #ifdef WITH_SMESH_CANCEL_COMPUTE
107   _compute_canceled = false;
108 #endif
109 }
110
111 //=============================================================================
112 /*!
113  *  
114  */
115 //=============================================================================
116
117 HexoticPlugin_Hexotic::~HexoticPlugin_Hexotic()
118 {
119   MESSAGE("HexoticPlugin_Hexotic::~HexoticPlugin_Hexotic");
120 }
121
122
123 #ifdef WITH_BLSURFPLUGIN
124 bool HexoticPlugin_Hexotic::CheckBLSURFHypothesis( SMESH_Mesh&         aMesh,
125                                                    const TopoDS_Shape& aShape )
126 {
127   // MESSAGE("HexoticPlugin_Hexotic::CheckBLSURFHypothesis");
128   _blsurfHypo = NULL;
129
130   std::list<const SMESHDS_Hypothesis*>::const_iterator itl;
131   const SMESHDS_Hypothesis* theHyp;
132
133   // If a BLSURF hypothesis is applied, get it
134   SMESH_HypoFilter blsurfFilter;
135   blsurfFilter.Init( blsurfFilter.HasName( "BLSURF_Parameters" ));
136   std::list<const SMESHDS_Hypothesis *> appliedHyps;
137   aMesh.GetHypotheses( aShape, blsurfFilter, appliedHyps, false );
138
139   if ( appliedHyps.size() > 0 ) {
140     itl = appliedHyps.begin();
141     theHyp = (*itl); // use only the first hypothesis
142     std::string hypName = theHyp->GetName();
143     if (hypName == "BLSURF_Parameters") {
144       _blsurfHypo = static_cast<const BLSURFPlugin_Hypothesis*> (theHyp);
145       ASSERT(_blsurfHypo);
146       return true;
147     }
148   }
149   return false;
150 }
151 #endif
152
153 //=============================================================================
154 /*!
155  *  
156  */
157 //=============================================================================
158
159 bool HexoticPlugin_Hexotic::CheckHypothesis( SMESH_Mesh&                          aMesh,
160                                              const TopoDS_Shape&                  aShape,
161                                              SMESH_Hypothesis::Hypothesis_Status& aStatus )
162 {
163   // MESSAGE("HexoticPlugin_Hexotic::CheckHypothesis");
164   _hypothesis = NULL;
165
166   std::list<const SMESHDS_Hypothesis*>::const_iterator itl;
167   const SMESHDS_Hypothesis* theHyp;
168
169   const std::list<const SMESHDS_Hypothesis*>& hyps = GetUsedHypothesis(aMesh, aShape, false);
170   int nbHyp = hyps.size();
171   if (!nbHyp) {
172     aStatus = SMESH_Hypothesis::HYP_OK;
173     return true;  // can work with no hypothesis
174   }
175
176   itl = hyps.begin();
177   theHyp = (*itl); // use only the first hypothesis
178
179   std::string hypName = theHyp->GetName();
180   if (hypName == "Hexotic_Parameters") {
181     _hypothesis = static_cast<const HexoticPlugin_Hypothesis*> (theHyp);
182     ASSERT(_hypothesis);
183     aStatus = SMESH_Hypothesis::HYP_OK;
184   }
185   else
186     aStatus = SMESH_Hypothesis::HYP_INCOMPATIBLE;
187   
188 #ifdef WITH_BLSURFPLUGIN
189   CheckBLSURFHypothesis(aMesh, aShape);
190 #endif
191   
192   return aStatus == SMESH_Hypothesis::HYP_OK;
193 }
194
195 //=======================================================================
196 //function : findShape
197 //purpose  :
198 //=======================================================================
199
200 static TopoDS_Shape findShape(SMDS_MeshNode**     t_Node,
201                               TopoDS_Shape        aShape,
202                               const TopoDS_Shape* t_Shape,
203                               double**            t_Box,
204                               const int           nShape)
205 {
206   double pntCoor[3];
207   int iShape, nbNode = 8;
208
209   for ( int i=0; i<3; i++ ) {
210     pntCoor[i] = 0;
211     for ( int j=0; j<nbNode; j++ ) {
212       if ( i == 0) pntCoor[i] += t_Node[j]->X();
213       if ( i == 1) pntCoor[i] += t_Node[j]->Y();
214       if ( i == 2) pntCoor[i] += t_Node[j]->Z();
215     }
216     pntCoor[i] /= nbNode;
217   }
218   gp_Pnt aPnt(pntCoor[0], pntCoor[1], pntCoor[2]);
219
220   if ( aShape.IsNull() ) aShape = t_Shape[0];
221   BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
222   if ( !(SC.State() == TopAbs_IN) ) {
223     aShape.Nullify();
224     for (iShape = 0; iShape < nShape && aShape.IsNull(); iShape++) {
225       if ( !( pntCoor[0] < t_Box[iShape][0] || t_Box[iShape][1] < pntCoor[0] ||
226               pntCoor[1] < t_Box[iShape][2] || t_Box[iShape][3] < pntCoor[1] ||
227               pntCoor[2] < t_Box[iShape][4] || t_Box[iShape][5] < pntCoor[2]) ) {
228         BRepClass3d_SolidClassifier SC (t_Shape[iShape], aPnt, Precision::Confusion());
229         if (SC.State() == TopAbs_IN)
230           aShape = t_Shape[iShape];
231       }
232     }
233   }
234   return aShape;
235 }
236
237 //=======================================================================
238 //function : findEdge
239 //purpose  :
240 //=======================================================================
241
242 static int findEdge(const SMDS_MeshNode* aNode,
243                     const SMESHDS_Mesh*  theMesh,
244                     const int            nEdge,
245                     const TopoDS_Shape*  t_Edge) {
246
247   TopoDS_Shape aPntShape, foundEdge;
248   TopoDS_Vertex aVertex;
249   gp_Pnt aPnt( aNode->X(), aNode->Y(), aNode->Z() );
250
251   int foundInd, ind;
252   double nearest = RealLast(), *t_Dist;
253   double epsilon = Precision::Confusion();
254
255   t_Dist = new double[ nEdge ];
256   aPntShape = BRepBuilderAPI_MakeVertex( aPnt ).Shape();
257   aVertex   = TopoDS::Vertex( aPntShape );
258
259   for ( ind=0; ind < nEdge; ind++ ) {
260     BRepExtrema_DistShapeShape aDistance ( aVertex, t_Edge[ind] );
261     t_Dist[ind] = aDistance.Value();
262     if ( t_Dist[ind] < nearest ) {
263       nearest   = t_Dist[ind];
264       foundEdge = t_Edge[ind];
265       foundInd  = ind;
266       if ( nearest < epsilon )
267         ind = nEdge;
268     }
269   }
270
271   delete [] t_Dist;
272   return theMesh->ShapeToIndex( foundEdge );
273 }
274
275 //=======================================================================
276 //function : getNbShape
277 //purpose  :
278 //=======================================================================
279
280 static int getNbShape(std::string aFile, std::string aString, int defaultValue=0) {
281   int number = defaultValue;
282   std::string aLine;
283   std::ifstream file(aFile.c_str());
284   while ( !file.eof() ) {
285     getline( file, aLine);
286     if ( aLine == aString ) {
287       getline( file, aLine);
288       std::istringstream stringFlux( aLine );
289       stringFlux >> number;
290       number = ( number + defaultValue + std::abs(number - defaultValue) ) / 2;
291       break;
292     }
293   }
294   file.close();
295   return number;
296 }
297
298 //=======================================================================
299 //function : countShape
300 //purpose  :
301 //=======================================================================
302
303 template < class Mesh, class Shape >
304 static int countShape( Mesh* mesh, Shape shape ) {
305   TopExp_Explorer expShape ( mesh->ShapeToMesh(), shape );
306   TopTools_MapOfShape mapShape;
307   int nbShape = 0;
308   for ( ; expShape.More(); expShape.Next() ) {
309     if (mapShape.Add(expShape.Current())) {
310       nbShape++;
311     }
312   }
313   return nbShape;
314 }
315
316 //=======================================================================
317 //function : getShape
318 //purpose  :
319 //=======================================================================
320
321 template < class Mesh, class Shape, class Tab >
322 void getShape(Mesh* mesh, Shape shape, Tab *t_Shape) {
323   TopExp_Explorer expShape ( mesh->ShapeToMesh(), shape );
324   TopTools_MapOfShape mapShape;
325   for ( int i=0; expShape.More(); expShape.Next() ) {
326     if (mapShape.Add(expShape.Current())) {
327       t_Shape[i] = expShape.Current();
328       i++;
329     }
330   }
331   return;
332 }
333
334 //=======================================================================
335 //function : printWarning
336 //purpose  :
337 //=======================================================================
338
339 static void printWarning(const int nbExpected, std::string aString, const int nbFound) {
340   cout << std::endl;
341   cout << "WARNING : " << nbExpected << " " << aString << " expected, Hexotic has found " << nbFound << std::endl;
342   cout << "=======" << std::endl;
343   cout << std::endl;
344   return;
345 }
346
347 //=======================================================================
348 //function : removeHexoticFiles
349 //purpose  :
350 //=======================================================================
351
352 static void removeHexoticFiles(TCollection_AsciiString file_In, TCollection_AsciiString file_Out) {
353   removeFile( file_In );
354   removeFile( file_Out );
355 }
356
357 //=======================================================================
358 //function : writeHexoticFile
359 //purpose  : 
360 //=======================================================================
361
362 static bool writeHexoticFile (std::ofstream&                       theFile,
363                               const SMESHDS_Mesh*                  theMesh,
364                               std::map <int,int>&                  theSmdsToHexoticIdMap,
365                               std::map <int,const SMDS_MeshNode*>& /*theHexoticIdToNodeMap*/,
366                               const TCollection_AsciiString&       Hexotic_In) {
367   cout << std::endl;
368   cout << "Creating Hexotic processed mesh file : " << Hexotic_In << std::endl;
369
370   int nbShape = 0;
371
372   TopExp_Explorer expface(theMesh->ShapeToMesh(), TopAbs_FACE);
373   for ( ; expface.More(); expface.Next() )
374     nbShape++;
375
376   int *tabID;
377   int *tabNodeId;
378   TopoDS_Shape *tabShape, aShape;
379
380   int shapeID;
381   bool idFound;
382
383   int nbVertices    = 0;
384   int nbTriangles   = 0;
385   const char* space = "  ";
386   int dummy_0D      = 0;
387   int dummy_2D;
388   int nbNodes = 0;
389
390   int aSmdsNodeID = 1;
391   const SMDS_MeshNode* aNode;
392   SMDS_NodeIteratorPtr itOnNode;
393
394   std::list< const SMDS_MeshElement* > faces;
395   std::list< const SMDS_MeshElement* >::iterator itListFace;
396   const SMDS_MeshElement* aFace;
397   SMESHDS_SubMesh* theSubMesh;
398   std::map<int,int>::const_iterator itOnSmdsNode;
399   SMDS_ElemIteratorPtr itOnSubNode, itOnSubFace;
400
401 // Writing SMESH points into Hexotic File
402
403   nbVertices = theMesh->NbNodes();
404
405   theFile << "MeshVersionFormatted 2" << std::endl;
406   theFile << std::endl;
407   theFile << "Dimension" << std::endl;
408   theFile << 3 << std::endl;
409   theFile << "# Set of mesh vertices" << std::endl;
410   theFile << "Vertices" << std::endl;
411   theFile << nbVertices << std::endl;
412
413   tabID     = new int[nbShape];
414   tabNodeId = new int[ nbVertices ];
415   tabShape  = new TopoDS_Shape[nbShape];
416
417   itOnNode = theMesh->nodesIterator();
418   while ( itOnNode->more() ) {
419       aNode = itOnNode->next();
420       dummy_0D = aNode->getshapeId();
421       tabNodeId[ aSmdsNodeID - 1 ] = 0;
422       idFound  = false;
423       for ( int j=0; j< aSmdsNodeID; j++ ) {
424         if ( dummy_0D == tabNodeId[j] ) {
425           idFound = true;
426           break;
427         }
428       }
429       if ( ! idFound )
430         tabNodeId[ aSmdsNodeID - 1 ] = dummy_0D;
431       theSmdsToHexoticIdMap.insert(std::map <int,int>::value_type( aNode->GetID(), aSmdsNodeID ));
432       aSmdsNodeID++;
433       theFile << aNode->X() << space << aNode->Y() << space << aNode->Z() << space << dummy_0D << std::endl;
434       }
435
436 // Writing SMESH faces into Hexotic File
437
438   ostringstream triaOut;
439   nbTriangles = 0;
440
441   expface.ReInit();
442   for ( int i = 0; expface.More(); expface.Next(), i++ ) {
443     tabID[i] = 0;
444     aShape   = expface.Current();
445     shapeID  = theMesh->ShapeToIndex( aShape );
446     idFound  = false;
447     for ( int j=0; j<=i; j++) {
448       if ( shapeID == tabID[j] ) {
449         idFound = true;
450         break;
451       }
452     }
453     if ( ! idFound ) {
454       tabID[i]    = shapeID;
455       tabShape[i] = aShape;
456     }
457   }
458   for ( int i=0; i<nbShape; i++ ) {
459     if ( ! (tabID[i] == 0) ) {
460       aShape      = tabShape[i];
461       shapeID     = tabID[i];
462       theSubMesh  = theMesh->MeshElements( aShape );
463       itOnSubFace = theSubMesh->GetElements();
464       while ( itOnSubFace->more() ) {
465         aFace    = itOnSubFace->next();
466         dummy_2D = shapeID;
467         
468         nbNodes = aFace->IsQuadratic() ? aFace->NbNodes()/2 : aFace->NbNodes();
469         if ( nbNodes == 3 ) // triangle
470         {
471           nbTriangles++;
472           for ( int i = 0; i < nbNodes; ++i )
473           {
474             aSmdsNodeID = aFace->GetNode( i )->GetID();
475             itOnSmdsNode = theSmdsToHexoticIdMap.find( aSmdsNodeID );
476             ASSERT( itOnSmdsNode != theSmdsToHexoticIdMap.end() );
477             triaOut << (*itOnSmdsNode).second << space;
478           }
479           triaOut << dummy_2D << std::endl;
480         }
481         else // polygon
482         {
483           int nbTria = nbNodes - 2, n0 = theSmdsToHexoticIdMap[ aFace->GetNode(0)->GetID() ];
484           nbTriangles += nbTria;
485           for ( int i = 0; i < nbTria; ++i )
486           {
487             triaOut << n0 << space;
488             triaOut << theSmdsToHexoticIdMap[ aFace->GetNode(i+1)->GetID() ] << space;
489             triaOut << theSmdsToHexoticIdMap[ aFace->GetNode(i+2)->GetID() ] << space;
490             triaOut << dummy_2D << std::endl;
491           }
492         }
493       }
494     }
495   }
496
497   theFile << std::endl;
498   theFile << "# Set of mesh triangles (v1,v2,v3,tag)" << std::endl;
499   theFile << "Triangles" << std::endl;
500   theFile << nbTriangles << std::endl;
501   theFile << triaOut.str() << std::endl;
502
503   theFile << std::endl;
504   theFile << "End" << std::endl;
505
506   cout << "Processed mesh file created, it contains :" << std::endl;
507   cout << "    " << nbVertices  << " vertices"  << std::endl;
508   cout << "    " << nbTriangles << " triangles" << std::endl;
509   cout << std::endl;
510
511   delete [] tabID;
512   delete [] tabNodeId;
513   delete [] tabShape;
514
515   return true;
516 }
517
518 //=======================================================================
519 //function : writeHexoticFile
520 //purpose  : 
521 //=======================================================================
522
523 static bool writeHexoticFile (std::ofstream&                       theFile,
524                               const SMESH_MesherHelper*            theHelper,
525                               std::map <int,int>&                  theSmdsToHexoticIdMap,
526                               const TCollection_AsciiString&       Hexotic_In)
527 {
528   cout << std::endl;
529   cout << "Creating Hexotic processed mesh file : " << Hexotic_In << std::endl;
530
531   int nbVertices    = 0;
532   int nbTriangles   = 0;
533   const char* space = "  ";
534   int dummy_0D      = 0;
535   int dummy_2D      = 0;
536
537   int aSmdsNodeID = 1;
538   const SMDS_MeshNode* aNode;
539   SMDS_NodeIteratorPtr itOnNode;
540
541   const SMDS_MeshElement* aFace;
542   std::map<int,int>::const_iterator itOnSmdsNode;
543   SMDS_ElemIteratorPtr itOnSubNode, itOnSubFace;
544
545   // Writing SMESH points into Hexotic File
546
547   nbVertices = theHelper->GetMeshDS()->NbNodes();
548
549   theFile << "MeshVersionFormatted 2" << std::endl;
550   theFile << std::endl;
551   theFile << "Dimension" << std::endl;
552   theFile << 3 << std::endl;
553   theFile << "# Set of mesh vertices" << std::endl;
554   theFile << "Vertices" << std::endl;
555   theFile << nbVertices << std::endl;
556
557   itOnNode = theHelper->GetMeshDS()->nodesIterator();
558   while ( itOnNode->more() )
559   {
560     aNode = itOnNode->next();
561     theSmdsToHexoticIdMap.insert(make_pair( aNode->GetID(), aSmdsNodeID ));
562     aSmdsNodeID++;
563     theFile << aNode->X() << space << aNode->Y() << space << aNode->Z() << space << dummy_0D << std::endl;
564   }
565
566   // Writing SMESH faces into Hexotic File
567
568   ostringstream triaOut;
569
570   itOnSubFace = theHelper->GetMeshDS()->elementsIterator(SMDSAbs_Face);
571   while ( itOnSubFace->more() )
572   {
573     aFace = itOnSubFace->next();
574     int nbNodes = aFace->IsQuadratic() ? aFace->NbNodes()/2 : aFace->NbNodes();
575     if ( nbNodes == 3 ) // triangle
576     {
577       nbTriangles++;
578       for ( int i = 0; i < nbNodes; ++i )
579       {
580         aSmdsNodeID = aFace->GetNode( i )->GetID();
581         itOnSmdsNode = theSmdsToHexoticIdMap.find( aSmdsNodeID );
582         ASSERT( itOnSmdsNode != theSmdsToHexoticIdMap.end() );
583         triaOut << (*itOnSmdsNode).second << space;
584       }
585       triaOut << dummy_2D << std::endl;
586     }
587     else // polygon
588     {
589       int nbTria = nbNodes - 2, n0 = theSmdsToHexoticIdMap[ aFace->GetNode(0)->GetID() ];
590       nbTriangles += nbTria;
591       for ( int i = 0; i < nbTria; ++i )
592       {
593         triaOut << n0 << space;
594         triaOut << theSmdsToHexoticIdMap[ aFace->GetNode(i+1)->GetID() ] << space;
595         triaOut << theSmdsToHexoticIdMap[ aFace->GetNode(i+2)->GetID() ] << space;
596         triaOut << dummy_2D << std::endl;
597       }
598     }
599   }
600
601   theFile << std::endl;
602   theFile << "# Set of mesh triangles (v1,v2,v3,tag)" << std::endl;
603   theFile << "Triangles" << std::endl;
604   theFile << nbTriangles << std::endl;
605   theFile << triaOut.str() << std::endl;
606
607   theFile << std::endl;
608   theFile << "End" << std::endl;
609
610   cout << "Processed mesh file created, it contains :" << std::endl;
611   cout << "    " << nbVertices  << " vertices"  << std::endl;
612   cout << "    " << nbTriangles << " triangles" << std::endl;
613   cout << std::endl;
614
615   return true;
616 }
617
618 //=======================================================================
619 //function : readResult
620 //purpose  : 
621 //=======================================================================
622
623 static bool readResult(std::string         theFile,
624 #ifdef WITH_SMESH_CANCEL_COMPUTE
625                        HexoticPlugin_Hexotic*  theAlgo,
626 #endif
627                        SMESHDS_Mesh*       theMesh,
628                        const int           nbShape,
629                        const TopoDS_Shape* tabShape,
630                        double**            tabBox)
631 {
632   // ---------------------------------
633   // Optimisation of the plugin ...
634   // Retrieve the correspondance edge --> shape
635   // (which is very costly) only when user
636   // has defined at least one group of edges
637   // which should be rare for a 3d mesh !
638   // ---------------------------------
639   
640   bool retrieve_edges = false;
641   const std::set<SMESHDS_GroupBase*>& aGroups = theMesh->GetGroups();
642   set<SMESHDS_GroupBase*>::const_iterator GrIt = aGroups.begin();
643   for (; GrIt != aGroups.end(); GrIt++)
644     {
645       SMESHDS_GroupBase* aGrp = *GrIt;
646       if ( !aGrp )
647         continue;
648       if ( aGrp->GetType() == SMDSAbs_Edge )
649         {
650           retrieve_edges = true;
651           break;
652         }
653     }
654   
655   // ---------------------------------
656   // Read generated elements and nodes
657   // ---------------------------------
658
659   TopoDS_Shape aShape;
660   TopoDS_Vertex aVertex;
661   std::string token;
662   int EndOfFile = 0, nbElem = 0, nField = 9, nbRef = 0;
663   int aHexoticNodeID = 0, shapeID, hexoticShapeID;
664   const int IdShapeRef = 2;
665   int *tabID, *tabRef, *nodeAssigne;
666   bool *tabDummy, hasDummy = false;
667   double epsilon = Precision::Confusion();
668   std::map <std::string,int> mapField;
669   SMDS_MeshNode** HexoticNode;
670   TopoDS_Shape *tabCorner, *tabEdge;
671
672   const int nbDomains = countShape( theMesh, TopAbs_SHELL );
673   const int holeID = -1;
674
675   // tabID    = new int[nbShape];
676   tabID    = new int[nbDomains];
677   tabRef   = new int[nField];
678   tabDummy = new bool[nField];
679
680   for (int i=0; i<nbDomains; i++)
681     tabID[i] = 0;
682   if ( nbDomains == 1 )
683     tabID[0] = theMesh->ShapeToIndex( tabShape[0] );
684
685   mapField["MeshVersionFormatted"] = 0; tabRef[0] = 0; tabDummy[0] = false;
686   mapField["Dimension"]            = 1; tabRef[1] = 0; tabDummy[1] = false;
687   mapField["Vertices"]             = 2; tabRef[2] = 3; tabDummy[2] = true;
688   mapField["Corners"]              = 3; tabRef[3] = 1; tabDummy[3] = false;
689   mapField["Edges"]                = 4; tabRef[4] = 2; tabDummy[4] = true;
690   mapField["Ridges"]               = 5; tabRef[5] = 1; tabDummy[5] = false;
691   mapField["Quadrilaterals"]       = 6; tabRef[6] = 4; tabDummy[6] = true;
692   mapField["Hexahedra"]            = 7; tabRef[7] = 8; tabDummy[7] = true;
693   mapField["End"]                  = 8; tabRef[8] = 0; tabDummy[0] = false;
694
695   SMDS_NodeIteratorPtr itOnHexoticInputNode = theMesh->nodesIterator();
696   while ( itOnHexoticInputNode->more() )
697     theMesh->RemoveNode( itOnHexoticInputNode->next() );
698
699   int nbVertices   = getNbShape(theFile, "Vertices");
700   int nbCorners    = getNbShape(theFile, "Corners", countShape( theMesh, TopAbs_VERTEX ));
701   int nbShapeEdge  = countShape( theMesh, TopAbs_EDGE );
702
703   tabCorner   = new TopoDS_Shape[ nbCorners ];
704   tabEdge     = new TopoDS_Shape[ nbShapeEdge ];
705   nodeAssigne = new int[ nbVertices + 1 ];
706   HexoticNode = new SMDS_MeshNode*[ nbVertices + 1 ];
707
708   getShape(theMesh, TopAbs_VERTEX, tabCorner);
709   getShape(theMesh, TopAbs_EDGE,   tabEdge);
710
711   MESSAGE("Read " << theFile << " file");
712   std::ifstream fileRes(theFile.c_str());
713   ASSERT(fileRes);
714
715   while ( EndOfFile == 0  ) {
716     int dummy;
717     fileRes >> token;
718
719     if (mapField.count(token)) {
720       nField   = mapField[token];
721       nbRef    = tabRef[nField];
722       hasDummy = tabDummy[nField];
723     }
724     else {
725       nField = -1;
726       nbRef = 0;
727     }
728
729     nbElem = 0;
730     if ( nField < (mapField.size() - 1) && nField >= 0 )
731       fileRes >> nbElem;
732
733     switch (nField) {
734       case 0: { // "MeshVersionFormatted"
735         MESSAGE(token << " " << nbElem);
736         break;
737       }
738       case 1: { // "Dimension"
739         MESSAGE("Mesh dimension " << nbElem << "D");
740         break;
741       }
742       case 2: { // "Vertices"
743         MESSAGE("Read " << nbElem << " " << token);
744         int aHexoticID;
745         double *coord;
746         SMDS_MeshNode * aHexoticNode;
747
748         coord = new double[nbRef];
749         for ( int iElem = 0; iElem < nbElem; iElem++ ) {
750 #ifdef WITH_SMESH_CANCEL_COMPUTE
751           if(theAlgo->computeCanceled())
752             {
753               return false;
754             }
755 #endif
756           aHexoticID = iElem + 1;
757           for ( int iCoord = 0; iCoord < 3; iCoord++ )
758             fileRes >> coord[ iCoord ];
759           fileRes >> dummy;
760           aHexoticNode = theMesh->AddNode(coord[0], coord[1], coord[2]);
761           HexoticNode[ aHexoticID ] = aHexoticNode;
762           nodeAssigne[ aHexoticID ] = 0;
763         }
764         delete [] coord;
765         break;
766       }
767       case 3: // "Corners"
768       case 4: // "Edges"
769       case 5: // "Ridges"
770       case 6: // "Quadrilaterals"
771       case 7: { // "Hexahedra"
772         MESSAGE("Read " << nbElem << " " << token);
773         SMDS_MeshNode** node;
774         int nodeDim, *nodeID;
775         SMDS_MeshElement * aHexoticElement = 0;
776
777         node   = new SMDS_MeshNode*[ nbRef ];
778         nodeID = new int[ nbRef ];
779         for ( int iElem = 0; iElem < nbElem; iElem++ ) {
780 #ifdef WITH_SMESH_CANCEL_COMPUTE
781           if(theAlgo->computeCanceled())
782             {
783               return false;
784             }
785 #endif
786           for ( int iRef = 0; iRef < nbRef; iRef++ ) {
787             fileRes >> aHexoticNodeID;                          // read nbRef aHexoticNodeID
788             node[ iRef ]   = HexoticNode[ aHexoticNodeID ];
789             nodeID[ iRef ] = aHexoticNodeID;
790           }
791           if ( hasDummy )
792             fileRes >> dummy;
793           switch (nField) {
794             case 3: { // "Corners"
795               nodeDim = 1;
796               gp_Pnt HexoticPnt ( node[0]->X(), node[0]->Y(), node[0]->Z() );
797               for ( int i=0; i<nbElem; i++ ) {
798                 aVertex = TopoDS::Vertex( tabCorner[i] );
799                 gp_Pnt aPnt = BRep_Tool::Pnt( aVertex );
800                 if ( aPnt.Distance( HexoticPnt ) < epsilon )
801                   break;
802               }
803               break;
804             }
805             case 4: { // "Edges"
806               nodeDim = 2;
807               aHexoticElement = theMesh->AddEdge( node[0], node[1] );
808               int iNode = 1;
809               if ( nodeAssigne[ nodeID[0] ] == 0 || nodeAssigne[ nodeID[0] ] == 2 )
810                 iNode = 0;
811               if(retrieve_edges)
812                 shapeID = findEdge( node[iNode], theMesh, nbShapeEdge, tabEdge );
813               else
814                 shapeID = 0;
815               break;
816             }
817             case 5: { // "Ridges"
818               break;
819             }
820             case 6: { // "Quadrilaterals"
821               nodeDim = 3;
822               aHexoticElement = theMesh->AddFace( node[0], node[1], node[2], node[3] );
823               shapeID = dummy;
824               break;
825             }
826             case 7: { // "Hexahedra"
827               nodeDim = 4;
828               if ( nbDomains > 1 ) {
829                 hexoticShapeID = dummy - IdShapeRef;
830                 if ( tabID[ hexoticShapeID ] == 0 ) {
831                   aShape = findShape(node, aShape, tabShape, tabBox, nbShape);
832                   shapeID = aShape.IsNull() ? holeID : theMesh->ShapeToIndex( aShape );
833                   tabID[ hexoticShapeID ] = shapeID;
834                 }
835                 else
836                   shapeID = tabID[ hexoticShapeID ];
837                 if ( iElem == (nbElem - 1) ) {
838                   int shapeAssociated = 0;
839                   for ( int i=0; i<nbDomains; i++ ) {
840                     if (tabID[i] > 0 )
841                       shapeAssociated += 1;
842                   }
843                   if ( shapeAssociated != nbShape )
844                     printWarning(nbShape, "domains", shapeAssociated);
845                 }
846               }
847               else {
848                 shapeID = tabID[0];
849               }
850               if ( shapeID != holeID )
851                 aHexoticElement = theMesh->AddVolume( node[0], node[3], node[2], node[1], node[4], node[7], node[6], node[5] );
852               break;
853             }
854           } // switch (nField)
855
856           if ( token != "Ridges" && ( shapeID > 0 || token == "Corners")) {
857             for ( int i=0; i<nbRef; i++ ) {
858               if ( nodeAssigne[ nodeID[i] ] == 0 ) {
859                 if      ( token == "Corners" )        theMesh->SetNodeOnVertex( node[0], aVertex );
860                 else if ( token == "Edges" )          theMesh->SetNodeOnEdge( node[i], shapeID );
861                 else if ( token == "Quadrilaterals" ) theMesh->SetNodeOnFace( node[i], shapeID );
862                 else if ( token == "Hexahedra" )      theMesh->SetNodeInVolume( node[i], shapeID );
863                 nodeAssigne[ nodeID[i] ] = nodeDim;
864               }
865             }
866             if ( token != "Corners" && aHexoticElement )
867               theMesh->SetMeshElementOnShape( aHexoticElement, shapeID );
868           }
869         }
870         delete [] node;
871         delete [] nodeID;
872         break;
873       }
874       case 8: { // "End"
875         EndOfFile = 1;
876         MESSAGE("End of " << theFile << " file");
877         break;
878       }
879       default: {
880         MESSAGE("Unknown Token: " << token);
881       }
882     }
883   }
884   cout << std::endl;
885
886   // remove nodes in holes
887   if ( nbDomains > 1 )
888   {
889     SMESHDS_SubMesh* subMesh;
890     for ( int i = 1; i <= nbVertices; ++i )
891       if ( HexoticNode[i]->NbInverseElements() == 0 )
892       {
893         subMesh =  HexoticNode[i]->getshapeId() > 0 ? theMesh->MeshElements(HexoticNode[i]->getshapeId() ) : 0;
894         theMesh->RemoveFreeNode( HexoticNode[i], subMesh, /*fromGroups=*/false );
895       }
896   }
897   delete [] tabID;
898   delete [] tabRef;
899   delete [] tabDummy;
900   delete [] tabCorner;
901   delete [] tabEdge;
902   delete [] nodeAssigne;
903   delete [] HexoticNode;
904   return true;
905 }
906
907
908 //=======================================================================
909 //function : readResult
910 //purpose  : 
911 //=======================================================================
912
913 static bool readResult(std::string theFile,
914 #ifdef WITH_SMESH_CANCEL_COMPUTE
915                        HexoticPlugin_Hexotic*  theAlgo,
916 #endif
917                        SMESH_MesherHelper* theHelper)
918 {
919   SMESHDS_Mesh* theMesh = theHelper->GetMeshDS();
920
921   // ---------------------------------
922   // Read generated elements and nodes
923   // ---------------------------------
924
925   std::string token;
926   const int nbField = 9;
927   int nField, EndOfFile = 0, nbElem = 0, nbRef = 0;
928   int aHexoticNodeID = 0, shapeID;
929   int tabRef[nbField], *nodeAssigne;
930   bool tabDummy[nbField], hasDummy = false;
931   std::map <std::string,int> mapField;
932   SMDS_MeshNode** HexoticNode;
933
934   mapField["MeshVersionFormatted"] = 0; tabRef[0] = 0; tabDummy[0] = false;
935   mapField["Dimension"]            = 1; tabRef[1] = 0; tabDummy[1] = false;
936   mapField["Vertices"]             = 2; tabRef[2] = 3; tabDummy[2] = true;
937   mapField["Corners"]              = 3; tabRef[3] = 1; tabDummy[3] = false;
938   mapField["Edges"]                = 4; tabRef[4] = 2; tabDummy[4] = true;
939   mapField["Ridges"]               = 5; tabRef[5] = 1; tabDummy[5] = false;
940   mapField["Quadrilaterals"]       = 6; tabRef[6] = 4; tabDummy[6] = true;
941   mapField["Hexahedra"]            = 7; tabRef[7] = 8; tabDummy[7] = true;
942   mapField["End"]                  = 8; tabRef[8] = 0; tabDummy[8] = false;
943
944   theHelper->GetMesh()->Clear();
945
946   int nbVertices = getNbShape(theFile, "Vertices");
947   HexoticNode = new SMDS_MeshNode*[ nbVertices + 1 ];
948   nodeAssigne = new int[ nbVertices + 1 ];
949
950   MESSAGE("Read " << theFile << " file");
951   std::ifstream fileRes(theFile.c_str());
952   ASSERT(fileRes);
953
954   while ( !EndOfFile  )
955   {
956     int dummy;
957     fileRes >> token;
958
959     if (mapField.count(token)) {
960       nField   = mapField[token];
961       nbRef    = tabRef[nField];
962       hasDummy = tabDummy[nField];
963     }
964     else {
965       nField = -1;
966       nbRef = 0;
967     }
968
969     nbElem = 0;
970     if ( nField < (mapField.size() - 1) && nField >= 0 )
971       fileRes >> nbElem;
972
973     switch (nField) {
974     case 0: { // "MeshVersionFormatted"
975       MESSAGE(token << " " << nbElem);
976       break;
977     }
978     case 1: { // "Dimension"
979       MESSAGE("Mesh dimension " << nbElem << "D");
980       break;
981     }
982     case 2: { // "Vertices"
983       MESSAGE("Read " << nbElem << " " << token);
984       int aHexoticID;
985       double coord[3];
986       SMDS_MeshNode * aHexoticNode;
987
988       for ( int iElem = 0; iElem < nbElem; iElem++ ) {
989 #ifdef WITH_SMESH_CANCEL_COMPUTE
990         if(theAlgo->computeCanceled())
991           {
992             return false;
993           }
994 #endif
995         aHexoticID = iElem + 1;
996         for ( int iCoord = 0; iCoord < 3; iCoord++ )
997           fileRes >> coord[ iCoord ];
998         fileRes >> dummy;
999         aHexoticNode = theMesh->AddNode(coord[0], coord[1], coord[2]);
1000         HexoticNode[ aHexoticID ] = aHexoticNode;
1001         nodeAssigne[ aHexoticID ] = 0;
1002       }
1003       break;
1004     }
1005     case 3: // "Corners"
1006     case 4: // "Edges"
1007     case 5: // "Ridges"
1008     case 6: // "Quadrilaterals"
1009     case 7: { // "Hexahedra"
1010       MESSAGE("Read " << nbElem << " " << token);
1011       std::vector< SMDS_MeshNode* > node( nbRef );
1012       std::vector< int >          nodeID( nbRef );
1013
1014       for ( int iElem = 0; iElem < nbElem; iElem++ )
1015       {
1016 #ifdef WITH_SMESH_CANCEL_COMPUTE
1017         if(theAlgo->computeCanceled())
1018           {
1019             return false;
1020           }
1021 #endif
1022         for ( int iRef = 0; iRef < nbRef; iRef++ )
1023         {
1024           fileRes >> aHexoticNodeID;                          // read nbRef aHexoticNodeID
1025           node  [ iRef ] = HexoticNode[ aHexoticNodeID ];
1026           nodeID[ iRef ] = aHexoticNodeID;
1027         }
1028         if ( hasDummy )
1029           fileRes >> dummy;
1030         switch (nField)
1031         {
1032         case 4: // "Edges"
1033           theHelper->AddEdge( node[0], node[1] ); break;
1034         case 6:  // "Quadrilaterals"
1035           theMesh->AddFace( node[0], node[1], node[2], node[3] ); break;
1036         case 7: // "Hexahedra"
1037           theHelper->AddVolume( node[0], node[3], node[2], node[1],
1038                                 node[4], node[7], node[6], node[5] ); break;
1039         default: continue;
1040         }
1041         if ( nField == 6 )
1042           for ( int iRef = 0; iRef < nbRef; iRef++ )
1043             nodeAssigne[ nodeID[ iRef ]] = 1;
1044       }
1045       break;
1046     }
1047     case 8: { // "End"
1048       EndOfFile = 1;
1049       MESSAGE("End of " << theFile << " file");
1050       break;
1051     }
1052     default: {
1053       MESSAGE("Unknown Token: " << token);
1054     }
1055     }
1056   }
1057   cout << std::endl;
1058
1059   shapeID = theHelper->GetSubShapeID();
1060   for ( int i = 0; i < nbVertices; ++i )
1061     if ( !nodeAssigne[ i+1 ])
1062       theMesh->SetNodeInVolume( HexoticNode[ i+1 ], shapeID );
1063
1064   delete [] HexoticNode;
1065   delete [] nodeAssigne;
1066   return true;
1067 }
1068
1069 //=============================================================================
1070 /*!
1071  * Pass parameters to Hexotic
1072  */
1073 //=============================================================================
1074
1075 void HexoticPlugin_Hexotic::SetParameters(const HexoticPlugin_Hypothesis* hyp) {
1076
1077   MESSAGE("HexoticPlugin_Hexotic::SetParameters");
1078   if (hyp) {
1079     _hexesMinLevel = hyp->GetHexesMinLevel();
1080     _hexesMaxLevel = hyp->GetHexesMaxLevel();
1081     _hexoticQuadrangles = hyp->GetHexoticQuadrangles();
1082     _hexoticIgnoreRidges = hyp->GetHexoticIgnoreRidges();
1083     _hexoticInvalidElements = hyp->GetHexoticInvalidElements();
1084     _hexoticSharpAngleThreshold = hyp->GetHexoticSharpAngleThreshold();
1085     _hexoticNbProc = hyp->GetHexoticNbProc();
1086     _hexoticWorkingDirectory = hyp->GetHexoticWorkingDirectory();
1087   }
1088   else {
1089     cout << std::endl;
1090     cout << "WARNING : The Hexotic default parameters are taken into account" << std::endl;
1091     cout << "=======" << std::endl;
1092     _hexesMinLevel = hyp->GetDefaultHexesMinLevel();
1093     _hexesMaxLevel = hyp->GetDefaultHexesMaxLevel();
1094     _hexoticQuadrangles = hyp->GetDefaultHexoticQuadrangles();
1095     _hexoticIgnoreRidges = hyp->GetDefaultHexoticIgnoreRidges();
1096     _hexoticInvalidElements = hyp->GetDefaultHexoticInvalidElements();
1097     _hexoticSharpAngleThreshold = hyp->GetDefaultHexoticSharpAngleThreshold();
1098     _hexoticNbProc = hyp->GetDefaultHexoticNbProc();
1099     _hexoticWorkingDirectory = hyp->GetDefaultHexoticWorkingDirectory();
1100   }
1101 }
1102
1103 //=======================================================================
1104 //function : getTmpDir
1105 //purpose  :
1106 //=======================================================================
1107
1108 static TCollection_AsciiString getTmpDir()
1109 {
1110   TCollection_AsciiString aTmpDir;
1111
1112   char *Tmp_dir = getenv("SALOME_TMP_DIR");
1113   if(Tmp_dir != NULL) {
1114     aTmpDir = Tmp_dir;
1115     #ifdef WIN32
1116     if(aTmpDir.Value(aTmpDir.Length()) != '\\') aTmpDir+='\\';
1117 #else
1118     if(aTmpDir.Value(aTmpDir.Length()) != '/') aTmpDir+='/';
1119 #endif
1120   }
1121   else {
1122 #ifdef WIN32
1123     aTmpDir = TCollection_AsciiString("C:\\");
1124 #else
1125     aTmpDir = TCollection_AsciiString("/tmp/");
1126 #endif
1127   }
1128   return aTmpDir;
1129 }
1130
1131 //================================================================================
1132 /*!
1133  * \brief Returns a command to run Hexotic mesher
1134  */
1135 //================================================================================
1136
1137 std::string HexoticPlugin_Hexotic::getHexoticCommand(const TCollection_AsciiString& Hexotic_In,
1138                                                      const TCollection_AsciiString& Hexotic_Out) const
1139 {
1140   cout << std::endl;
1141   cout << "Hexotic execution..." << std::endl;
1142   cout << _name << " parameters :" << std::endl;
1143   cout << "    " << _name << " Segments Min Level = " << _hexesMinLevel << std::endl;
1144   cout << "    " << _name << " Segments Max Level = " << _hexesMaxLevel << std::endl;
1145   cout << "    " << "Salome Quadrangles : " << (_hexoticQuadrangles ? "yes":"no") << std::endl;
1146   cout << "    " << "Hexotic can ignore ridges : " << (_hexoticIgnoreRidges ? "yes":"no") << std::endl;
1147   cout << "    " << "Hexotic authorize invalide elements : " << ( _hexoticInvalidElements ? "yes":"no") << std::endl;
1148   cout << "    " << _name << " Sharp angle threshold = " << _hexoticSharpAngleThreshold << " degrees" << std::endl;
1149   cout << "    " << _name << " Number of threads = " << _hexoticNbProc << std::endl;
1150   cout << "    " << _name << " Working directory = \"" << _hexoticWorkingDirectory << "\"" << std::endl;
1151
1152   TCollection_AsciiString run_Hexotic( "hexotic" );
1153
1154   TCollection_AsciiString minl = " -minl ", maxl = " -maxl ", angle = " -ra ";
1155   TCollection_AsciiString in   = " -in ",   out  = " -out ";
1156   TCollection_AsciiString ignoreRidges = " -nr ", invalideElements = " -inv ";
1157   TCollection_AsciiString subdom = " -sd ", sharp = " -sharp ";
1158   TCollection_AsciiString proc = " -nproc ";
1159
1160   TCollection_AsciiString minLevel, maxLevel, sharpAngle, mode, subdivision, nbproc;
1161   minLevel = _hexesMinLevel;
1162   maxLevel = _hexesMaxLevel;
1163   sharpAngle = _hexoticSharpAngleThreshold;
1164   mode = 4;
1165   subdivision = 3;
1166   nbproc = _hexoticNbProc;
1167
1168   if (_hexoticIgnoreRidges)
1169     run_Hexotic +=  ignoreRidges;
1170
1171   if (_hexoticInvalidElements)
1172     run_Hexotic +=  invalideElements;
1173
1174   run_Hexotic += angle + sharpAngle + minl + minLevel + maxl + maxLevel + in + Hexotic_In + out + Hexotic_Out;
1175   run_Hexotic += subdom + mode;
1176   run_Hexotic += proc + nbproc;
1177   //     run_Hexotic += subdom + mode + invalideElements;
1178   //     run_Hexotic += subdom + mode + ignoreRidges;
1179   //     run_Hexotic += subdom + mode + sharp + subdivision;
1180
1181   return run_Hexotic.ToCString();
1182 }
1183
1184 //=============================================================================
1185 /*!
1186  * Here we are going to use the Hexotic mesher
1187  */
1188 //=============================================================================
1189
1190 bool HexoticPlugin_Hexotic::Compute(SMESH_Mesh&          theMesh,
1191                                      const TopoDS_Shape& theShape)
1192 {
1193 #ifdef WITH_SMESH_CANCEL_COMPUTE
1194   _compute_canceled = false;
1195 #endif
1196   bool Ok = true;
1197   SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1198   TCollection_AsciiString hexahedraMessage;
1199
1200   if (_iShape == 0 && _nbShape == 0) {
1201     _nbShape = countShape( meshDS, TopAbs_SOLID );  // we count the number of shapes
1202     //_tabNode = new SMDS_MeshNode*[_nbShape];        // we declare the size of the node array
1203   }
1204
1205   // to prevent from displaying error message after computing,
1206   // SetIsAlwaysComputed( true ) to empty sub-meshes
1207   for ( int i = 0; i < _nbShape; ++i )
1208     if ( SMESH_subMesh* sm = theMesh.GetSubMeshContaining( theShape ))
1209     {
1210       SMESH_subMeshIteratorPtr smIt = sm->getDependsOnIterator(/*includeSelf=*/true,
1211                                                                /*complexShapeFirst=*/false);
1212       while ( smIt->more() )
1213       {
1214         sm = smIt->next();
1215         if ( !sm->IsMeshComputed() )
1216           sm->SetIsAlwaysComputed( true );
1217       }
1218     }
1219
1220   _iShape++;
1221
1222   if (_iShape == _nbShape ) {
1223
1224     // for (int i=0; i<_nbShape; i++)        // we destroy the (_nbShape - 1) nodes created and used
1225     //   meshDS->RemoveNode( _tabNode[i] );  // to simulate successful mesh computing.
1226     // delete [] _tabNode;
1227
1228     // create bounding box for each shape of the compound
1229
1230     int iShape = 0;
1231     TopoDS_Shape *tabShape;
1232     double **tabBox;
1233
1234     tabShape = new TopoDS_Shape[_nbShape];
1235     tabBox   = new double*[_nbShape];
1236     for (int i=0; i<_nbShape; i++)
1237       tabBox[i] = new double[6];
1238     double Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
1239
1240     TopExp_Explorer expBox (meshDS->ShapeToMesh(), TopAbs_SOLID);
1241     for (; expBox.More(); expBox.Next()) {
1242       tabShape[iShape] = expBox.Current();
1243       Bnd_Box BoundingBox;
1244       BRepBndLib::Add(expBox.Current(), BoundingBox);
1245       BoundingBox.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
1246       tabBox[iShape][0] = Xmin; tabBox[iShape][1] = Xmax;
1247       tabBox[iShape][2] = Ymin; tabBox[iShape][3] = Ymax;
1248       tabBox[iShape][4] = Zmin; tabBox[iShape][5] = Zmax;
1249       iShape++;
1250     }
1251
1252     SetParameters(_hypothesis);
1253
1254 //     TCollection_AsciiString aTmpDir = getTmpDir();
1255     TCollection_AsciiString aTmpDir = TCollection_AsciiString(_hexoticWorkingDirectory.c_str());
1256 #ifdef WIN32
1257     if ( aTmpDir.Value(aTmpDir.Length()) != '\\' ) aTmpDir += '\\';
1258 #else
1259     if ( aTmpDir.Value(aTmpDir.Length()) != '/' ) aTmpDir += '/';
1260 #endif
1261     TCollection_AsciiString Hexotic_In(""), Hexotic_Out;
1262     TCollection_AsciiString modeFile_In( "chmod 666 " ), modeFile_Out( "chmod 666 " );
1263     TCollection_AsciiString aLogFileName = aTmpDir + "Hexotic.log";    // log
1264
1265     std::map <int,int> aSmdsToHexoticIdMap;
1266     std::map <int,const SMDS_MeshNode*> aHexoticIdToNodeMap;
1267
1268     Hexotic_Out = aTmpDir + "Hexotic_Out.mesh";
1269 #ifdef WITH_BLSURFPLUGIN
1270     bool defaultInputFile = true;
1271     if (_blsurfHypo && !_blsurfHypo->GetQuadAllowed()) {
1272       Hexotic_In = TCollection_AsciiString(_blsurfHypo->GetGMFFile().c_str());
1273       if (Hexotic_In != "")
1274         defaultInputFile = false;
1275     }
1276     if (defaultInputFile) {
1277 #endif
1278       Hexotic_In  = aTmpDir + "Hexotic_In.mesh";
1279       removeHexoticFiles(Hexotic_In, Hexotic_Out);
1280       std::ofstream HexoticFile (Hexotic_In.ToCString(), std::ios::out);
1281       Ok = ( writeHexoticFile(HexoticFile, meshDS, aSmdsToHexoticIdMap, aHexoticIdToNodeMap, Hexotic_In) );
1282       HexoticFile.close();
1283 #ifdef WITH_BLSURFPLUGIN
1284     }
1285     else {
1286       removeFile( Hexotic_Out );
1287     }
1288 #endif
1289     
1290
1291     std::string run_Hexotic = getHexoticCommand(Hexotic_In, Hexotic_Out);
1292     run_Hexotic += std::string(" 1>") + aLogFileName.ToCString();  // dump into file
1293
1294     cout << std::endl;
1295     cout << "Hexotic command : " << run_Hexotic << std::endl;
1296
1297     
1298 //     removeHexoticFiles(Hexotic_In, Hexotic_Out);
1299 // 
1300 //     std::ofstream HexoticFile (Hexotic_In.ToCString(), std::ios::out);
1301 // 
1302 //     Ok = ( writeHexoticFile(HexoticFile, meshDS, aSmdsToHexoticIdMap, aHexoticIdToNodeMap, Hexotic_In) );
1303
1304 //     HexoticFile.close();
1305     modeFile_In += Hexotic_In;
1306     system( modeFile_In.ToCString() );
1307     aSmdsToHexoticIdMap.clear();
1308     aHexoticIdToNodeMap.clear();
1309
1310     MESSAGE("HexoticPlugin_Hexotic::Compute");
1311
1312     system( run_Hexotic.data() );
1313
1314     // --------------
1315     // read a result
1316     // --------------
1317
1318     std::ifstream fileRes( Hexotic_Out.ToCString() );
1319     modeFile_Out += Hexotic_Out;
1320     system( modeFile_Out.ToCString() );
1321     if ( ! fileRes.fail() ) {
1322       Ok = readResult( Hexotic_Out.ToCString(),
1323 #ifdef WITH_SMESH_CANCEL_COMPUTE
1324                        this,
1325 #endif
1326                        meshDS, _nbShape, tabShape, tabBox );
1327       if(Ok)
1328         hexahedraMessage = "success";
1329       else
1330         hexahedraMessage = "failed";
1331     }
1332     else {
1333       hexahedraMessage = "failed";
1334       cout << "Problem with Hexotic output file " << Hexotic_Out.ToCString() << std::endl;
1335       Ok = false;
1336       // analyse log file
1337       SMESH_File logFile( aLogFileName.ToCString() );
1338       if ( !logFile.eof() )
1339       {
1340         char msgLic[] = " Dlim ";
1341         const char* fileBeg = logFile.getPos(), *fileEnd = fileBeg + logFile.size();
1342         if ( std::search( fileBeg, fileEnd, msgLic, msgLic+strlen(msgLic)) != fileEnd )
1343           error("Licence problems.");
1344       }
1345     }
1346     cout << "Hexahedra meshing " << hexahedraMessage << std::endl;
1347     cout << std::endl;
1348
1349     delete [] tabShape;
1350     for (int i=0; i<_nbShape; i++)
1351       delete [] tabBox[i];
1352     delete [] tabBox;
1353     _nbShape = 0;
1354     _iShape  = 0;
1355   }
1356 #ifdef WITH_SMESH_CANCEL_COMPUTE
1357   if(_compute_canceled)
1358     return error(SMESH_Comment("interruption initiated by user"));
1359 #endif
1360   return Ok;
1361 }
1362
1363 //=============================================================================
1364 /*!
1365  * \brief Computes mesh without geometry
1366  *  \param aMesh - the mesh
1367  *  \param aHelper - helper that must be used for adding elements to \aaMesh
1368  *  \retval bool - is a success
1369  *
1370  * The method is called if ( !aMesh->HasShapeToMesh() )
1371  */
1372 //=============================================================================
1373
1374 bool HexoticPlugin_Hexotic::Compute(SMESH_Mesh & aMesh, SMESH_MesherHelper* aHelper)
1375 {
1376 #ifdef WITH_SMESH_CANCEL_COMPUTE
1377   _compute_canceled = false;
1378 #endif
1379   bool Ok = true;
1380   TCollection_AsciiString hexahedraMessage;
1381
1382   SetParameters(_hypothesis);
1383
1384   TCollection_AsciiString aTmpDir = getTmpDir();
1385   TCollection_AsciiString Hexotic_In, Hexotic_Out;
1386   TCollection_AsciiString modeFile_In( "chmod 666 " ), modeFile_Out( "chmod 666 " );
1387   TCollection_AsciiString aLogFileName = aTmpDir + "Hexotic.log";    // log
1388
1389   std::map <int,int> aSmdsToHexoticIdMap;
1390   std::map <int,const SMDS_MeshNode*> aHexoticIdToNodeMap;
1391
1392   Hexotic_In  = aTmpDir + "Hexotic_In.mesh";
1393   Hexotic_Out = aTmpDir + "Hexotic_Out.mesh";
1394
1395   std::string run_Hexotic = getHexoticCommand(Hexotic_In, Hexotic_Out);
1396   run_Hexotic += std::string(" 1>") + aLogFileName.ToCString();  // dump into file
1397
1398   cout << std::endl;
1399   cout << "Hexotic command : " << run_Hexotic << std::endl;
1400
1401   removeHexoticFiles(Hexotic_In, Hexotic_Out);
1402
1403   std::ofstream HexoticFile (Hexotic_In.ToCString(), std::ios::out);
1404
1405   Ok = ( writeHexoticFile(HexoticFile, aHelper, aSmdsToHexoticIdMap, Hexotic_In) );
1406
1407   HexoticFile.close();
1408   modeFile_In += Hexotic_In;
1409   system( modeFile_In.ToCString() );
1410   aSmdsToHexoticIdMap.clear();
1411   aHexoticIdToNodeMap.clear();
1412
1413   MESSAGE("HexoticPlugin_Hexotic::Compute");
1414
1415   system( run_Hexotic.data() );
1416
1417   // --------------
1418   // read a result
1419   // --------------
1420
1421   std::ifstream fileRes( Hexotic_Out.ToCString() );
1422   modeFile_Out += Hexotic_Out;
1423   system( modeFile_Out.ToCString() );
1424   if ( ! fileRes.fail() ) {
1425     Ok = readResult( Hexotic_Out.ToCString(),
1426 #ifdef WITH_SMESH_CANCEL_COMPUTE
1427                      this,
1428 #endif
1429                      aHelper );
1430     if(Ok)
1431       hexahedraMessage = "success";
1432     else
1433       hexahedraMessage = "failed";
1434   }
1435   else {
1436     hexahedraMessage = "failed";
1437     cout << "Problem with Hexotic output file " << Hexotic_Out << std::endl;
1438     // analyse log file
1439     SMESH_File logFile( aLogFileName.ToCString() );
1440     if ( !logFile.eof() )
1441     {
1442       char msgLic[] = " Dlim ";
1443       const char* fileBeg = logFile.getPos(), *fileEnd = fileBeg + logFile.size();
1444       if ( std::search( fileBeg, fileEnd, msgLic, msgLic+strlen(msgLic)) != fileEnd )
1445         return error("Licence problems.");
1446     }
1447     return error(SMESH_Comment("Problem with Hexotic output file ")<<Hexotic_Out);
1448   }
1449   cout << "Hexahedra meshing " << hexahedraMessage << std::endl;
1450   cout << std::endl;
1451
1452 #ifdef WITH_SMESH_CANCEL_COMPUTE
1453   if(_compute_canceled)
1454     return error(SMESH_Comment("interruption initiated by user"));
1455 #endif
1456   return Ok;
1457 }
1458
1459 //=============================================================================
1460 /*!
1461  *
1462  */
1463 //=============================================================================
1464
1465 bool HexoticPlugin_Hexotic::Evaluate(SMESH_Mesh& aMesh,
1466                                      const TopoDS_Shape& aShape,
1467                                      MapShapeNbElems& aResMap)
1468 {
1469   std::vector<int> aResVec(SMDSEntity_Last);
1470   for(int i=SMDSEntity_Node; i<SMDSEntity_Last; i++) aResVec[i] = 0;
1471   SMESH_subMesh * sm = aMesh.GetSubMesh(aShape);
1472   aResMap.insert(std::make_pair(sm,aResVec));
1473
1474   SMESH_ComputeErrorPtr& smError = sm->GetComputeError();
1475   smError.reset( new SMESH_ComputeError(COMPERR_ALGO_FAILED,"Evaluation is not implemented",this));
1476
1477   return true;
1478 }
1479
1480 #ifdef WITH_SMESH_CANCEL_COMPUTE
1481 void HexoticPlugin_Hexotic::CancelCompute()
1482 {
1483   _compute_canceled = true;
1484 #ifdef WNT
1485 #else
1486   TCollection_AsciiString aTmpDir = getTmpDir();
1487   TCollection_AsciiString Hexotic_In = aTmpDir + "Hexotic_In.mesh";
1488   TCollection_AsciiString cmd = TCollection_AsciiString("ps ux | grep ") + Hexotic_In;
1489   cmd += TCollection_AsciiString(" | grep -v grep | awk '{print $2}' | xargs kill -9 > /dev/null 2>&1");
1490   system( cmd.ToCString() );
1491 #endif
1492 }
1493 #endif