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