]> SALOME platform Git repositories - plugins/ghs3dplugin.git/blob - src/GHS3DPlugin_GHS3D.cxx
Salome HOME
PAL19680 Meshers: BLSURF, GHS3D and holed shapes
[plugins/ghs3dplugin.git] / src / GHS3DPlugin_GHS3D.cxx
1 // Copyright (C) 2005  CEA/DEN, EDF R&D, OPEN CASCADE, PRINCIPIA 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 // File      : GHS3DPlugin_GHS3D.cxx
21 // Created   : 
22 // Author    : Edward AGAPOV, modified by Lioka RAZAFINDRAZAKA (CEA) 09/02/2007
23 // Project   : SALOME
24 // Copyright : CEA 2003
25 // $Header$
26 //=============================================================================
27 using namespace std;
28
29 #include "GHS3DPlugin_GHS3D.hxx"
30 #include "GHS3DPlugin_Hypothesis.hxx"
31
32 #include "SMESH_Gen.hxx"
33 #include "SMESH_Mesh.hxx"
34 #include "SMESH_Comment.hxx"
35 #include "SMESH_MesherHelper.hxx"
36
37 #include "SMDS_MeshElement.hxx"
38 #include "SMDS_MeshNode.hxx"
39
40 #include <TopExp.hxx>
41 #include <TopExp_Explorer.hxx>
42 #include <OSD_File.hxx>
43
44 #include "utilities.h"
45
46 #ifndef WIN32
47 #include <sys/sysinfo.h>
48 #endif
49
50 //#include <Standard_Stream.hxx>
51
52 #include <BRepGProp.hxx>
53 #include <BRepBndLib.hxx>
54 #include <BRepClass_FaceClassifier.hxx>
55 #include <BRepClass3d_SolidClassifier.hxx>
56 #include <TopAbs.hxx>
57 #include <Bnd_Box.hxx>
58 #include <GProp_GProps.hxx>
59 #include <Precision.hxx>
60
61 #define castToNode(n) static_cast<const SMDS_MeshNode *>( n );
62
63 #ifdef _DEBUG_
64 #define DUMP(txt) \
65 //  cout << txt
66 #else
67 #define DUMP(txt)
68 #endif
69
70 extern "C"
71 {
72 #ifndef WNT
73 #include <unistd.h>
74 #include <sys/mman.h>
75 #endif
76 #include <sys/stat.h>
77 #include <fcntl.h>
78 }
79
80 //=============================================================================
81 /*!
82  *  
83  */
84 //=============================================================================
85
86 GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D(int hypId, int studyId, SMESH_Gen* gen)
87   : SMESH_3D_Algo(hypId, studyId, gen)
88 {
89   MESSAGE("GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D");
90   _name = "GHS3D_3D";
91   _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);// 1 bit /shape type
92   _onlyUnaryInput = false; // Compute() will be called on a compound of solids
93   _iShape=0;
94   _nbShape=0;
95   _compatibleHypothesis.push_back("GHS3D_Parameters");
96 }
97
98 //=============================================================================
99 /*!
100  *  
101  */
102 //=============================================================================
103
104 GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D()
105 {
106   MESSAGE("GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D");
107 }
108
109 //=============================================================================
110 /*!
111  *  
112  */
113 //=============================================================================
114
115 bool GHS3DPlugin_GHS3D::CheckHypothesis ( SMESH_Mesh&         aMesh,
116                                           const TopoDS_Shape& aShape,
117                                           Hypothesis_Status&  aStatus )
118 {
119   aStatus = SMESH_Hypothesis::HYP_OK;
120
121   // there is only one compatible Hypothesis so far
122   _hyp = 0;
123   _keepFiles = false;
124   const list <const SMESHDS_Hypothesis * >& hyps = GetUsedHypothesis(aMesh, aShape);
125   if ( !hyps.empty() )
126     _hyp = static_cast<const GHS3DPlugin_Hypothesis*> ( hyps.front() );
127   if ( _hyp )
128     _keepFiles = _hyp->GetKeepFiles();
129
130   return true;
131 }
132
133 //=======================================================================
134 //function : findShape
135 //purpose  : 
136 //=======================================================================
137
138 static TopoDS_Shape findShape(const SMDS_MeshNode *aNode[],
139                               TopoDS_Shape        aShape,
140                               const TopoDS_Shape  shape[],
141                               double**            box,
142                               const int           nShape,
143                               TopAbs_State *      state = 0) {
144   double *pntCoor;
145   int iShape, nbNode = 4;
146
147   pntCoor = new double[3];
148   for ( int i=0; i<3; i++ ) {
149     pntCoor[i] = 0;
150     for ( int j=0; j<nbNode; j++ ) {
151       if ( i == 0) pntCoor[i] += aNode[j]->X();
152       if ( i == 1) pntCoor[i] += aNode[j]->Y();
153       if ( i == 2) pntCoor[i] += aNode[j]->Z();
154     }
155     pntCoor[i] /= nbNode;
156   }
157
158   gp_Pnt aPnt(pntCoor[0], pntCoor[1], pntCoor[2]);
159   BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
160   if (state) *state = SC.State();
161   if ( SC.State() != TopAbs_IN || aShape.IsNull() || aShape.ShapeType() != TopAbs_SOLID) {
162     for (iShape = 0; iShape < nShape; iShape++) {
163       aShape = shape[iShape];
164       if ( !( pntCoor[0] < box[iShape][0] || box[iShape][1] < pntCoor[0] ||
165               pntCoor[1] < box[iShape][2] || box[iShape][3] < pntCoor[1] ||
166               pntCoor[2] < box[iShape][4] || box[iShape][5] < pntCoor[2]) ) {
167         BRepClass3d_SolidClassifier SC (aShape, aPnt, Precision::Confusion());
168         if (state) *state = SC.State();
169         if (SC.State() == TopAbs_IN)
170           break;
171       }
172     }
173   }
174   return aShape;
175 }
176
177 //=======================================================================
178 //function : readMapIntLine
179 //purpose  : 
180 //=======================================================================
181
182 static char* readMapIntLine(char* ptr, int tab[]) {
183   long int intVal;
184   cout << endl;
185
186   for ( int i=0; i<17; i++ ) {
187     intVal = strtol(ptr, &ptr, 10);
188     if ( i < 3 )
189       tab[i] = intVal;
190   }
191   return ptr;
192 }
193
194 //=======================================================================
195 //function : readLine
196 //purpose  : 
197 //=======================================================================
198
199 #define GHS3DPlugin_BUFLENGTH 256
200 #define GHS3DPlugin_ReadLine(aPtr,aBuf,aFile,aLineNb) \
201 {  aPtr = fgets( aBuf, GHS3DPlugin_BUFLENGTH - 2, aFile ); aLineNb++; DUMP(endl); }
202
203 //=======================================================================
204 //function : countShape
205 //purpose  :
206 //=======================================================================
207
208 template < class Mesh, class Shape >
209 static int countShape( Mesh* mesh, Shape shape ) {
210   TopExp_Explorer expShape ( mesh->ShapeToMesh(), shape );
211   int nbShape = 0;
212   for ( ; expShape.More(); expShape.Next() ) {
213       nbShape++;
214   }
215   return nbShape;
216 }
217
218 //=======================================================================
219 //function : writeFaces
220 //purpose  : 
221 //=======================================================================
222
223 static bool writeFaces (ofstream &            theFile,
224                         SMESHDS_Mesh *        theMesh,
225                         const map <int,int> & theSmdsToGhs3dIdMap)
226 {
227   // record structure:
228   //
229   // NB_ELEMS DUMMY_INT
230   // Loop from 1 to NB_ELEMS
231   // NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
232
233   int nbShape = countShape( theMesh, TopAbs_FACE );
234
235   int *tabID;             tabID    = new int[nbShape];
236   TopoDS_Shape *tabShape; tabShape = new TopoDS_Shape[nbShape];
237   TopoDS_Shape aShape;
238   SMESHDS_SubMesh* theSubMesh;
239   const SMDS_MeshElement* aFace;
240   const char* space    = "  ";
241   const int   dummyint = 0;
242   map<int,int>::const_iterator itOnMap;
243   SMDS_ElemIteratorPtr itOnSubMesh, itOnSubFace;
244   int shapeID, nbNodes, aSmdsID;
245   bool idFound;
246
247   cout << "    " << theMesh->NbFaces() << " shapes of 2D dimension" << endl;
248   cout << endl;
249
250   theFile << space << theMesh->NbFaces() << space << dummyint << endl;
251
252   TopExp_Explorer expface( theMesh->ShapeToMesh(), TopAbs_FACE );
253   for ( int i = 0; expface.More(); expface.Next(), i++ ) {
254     tabID[i] = 0;
255     aShape   = expface.Current();
256     shapeID  = theMesh->ShapeToIndex( aShape );
257     idFound  = false;
258     for ( int j=0; j<=i; j++) {
259       if ( shapeID == tabID[j] ) {
260         idFound = true;
261         break;
262       }
263     }
264     if ( ! idFound ) {
265       tabID[i]    = shapeID;
266       tabShape[i] = aShape;
267     }
268   }
269   for ( int i =0; i < nbShape; i++ ) {
270     if ( tabID[i] != 0 ) {
271       aShape      = tabShape[i];
272       shapeID     = tabID[i];
273       theSubMesh  = theMesh->MeshElements( aShape );
274       if ( !theSubMesh ) continue;
275       itOnSubMesh = theSubMesh->GetElements();
276       while ( itOnSubMesh->more() ) {
277         aFace   = itOnSubMesh->next();
278         nbNodes = aFace->NbNodes();
279
280         theFile << space << nbNodes;
281
282         itOnSubFace = aFace->nodesIterator();
283         while ( itOnSubFace->more() ) {
284           // find GHS3D ID
285           aSmdsID = itOnSubFace->next()->GetID();
286           itOnMap = theSmdsToGhs3dIdMap.find( aSmdsID );
287           ASSERT( itOnMap != theSmdsToGhs3dIdMap.end() );
288
289           theFile << space << (*itOnMap).second;
290         }
291
292         // (NB_NODES + 1) times: DUMMY_INT
293         for ( int j=0; j<=nbNodes; j++)
294           theFile << space << dummyint;
295
296         theFile << endl;
297       }
298     }
299   }
300
301   delete [] tabID;
302   delete [] tabShape;
303
304   return true;
305 }
306
307 //=======================================================================
308 //function : writeFaces
309 //purpose  : Write Faces in case if generate 3D mesh w/o geometry
310 //=======================================================================
311
312 static bool writeFaces (ofstream &            theFile,
313                         SMESHDS_Mesh *        theMesh,
314                         vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
315 {
316   // record structure:
317   //
318   // NB_ELEMS DUMMY_INT
319   // Loop from 1 to NB_ELEMS
320   //   NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
321
322
323   int nbFaces = 0;
324   list< const SMDS_MeshElement* > faces;
325   list< const SMDS_MeshElement* >::iterator f;
326   map< const SMDS_MeshNode*,int >::iterator it;
327   SMDS_ElemIteratorPtr nodeIt;
328   const SMDS_MeshElement* elem;
329   int nbNodes;
330
331   const char* space    = "  ";
332   const int   dummyint = 0;
333
334   //get all faces from mesh
335   SMDS_FaceIteratorPtr eIt = theMesh->facesIterator();
336   while ( eIt->more() ) {
337     const SMDS_MeshElement* elem = eIt->next();
338     if ( !elem )
339       return false;
340     faces.push_back( elem );
341     nbFaces++;
342   }
343
344   if ( nbFaces == 0 )
345     return false;
346   
347   cout << "The initial 2D mesh contains " << nbFaces << " faces and ";
348
349   // NB_ELEMS DUMMY_INT
350   theFile << space << nbFaces << space << dummyint << endl;
351
352   // Loop from 1 to NB_ELEMS
353
354   map<const SMDS_MeshNode*,int> aNodeToGhs3dIdMap;
355   f = faces.begin();
356   for ( ; f != faces.end(); ++f )
357   {
358     // NB_NODES PER FACE
359     elem = *f;
360     nbNodes = elem->NbNodes();
361     theFile << space << nbNodes;
362
363     // NODE_NB_1 NODE_NB_2 ...
364     nodeIt = elem->nodesIterator();
365     while ( nodeIt->more() )
366     {
367       // find GHS3D ID
368       const SMDS_MeshNode* node = castToNode( nodeIt->next() );
369       int newId = aNodeToGhs3dIdMap.size() + 1; // ghs3d ids count from 1
370       it = aNodeToGhs3dIdMap.insert( make_pair( node, newId )).first;
371       theFile << space << it->second;
372     }
373
374     // (NB_NODES + 1) times: DUMMY_INT
375     for ( int i=0; i<=nbNodes; i++)
376       theFile << space << dummyint;
377
378     theFile << endl;
379   }
380
381   // put nodes to theNodeByGhs3dId vector
382   theNodeByGhs3dId.resize( aNodeToGhs3dIdMap.size() );
383   map<const SMDS_MeshNode*,int>::const_iterator n2id = aNodeToGhs3dIdMap.begin();
384   for ( ; n2id != aNodeToGhs3dIdMap.end(); ++ n2id)
385   {
386     theNodeByGhs3dId[ n2id->second - 1 ] = n2id->first; // ghs3d ids count from 1
387   }
388
389   return true;
390 }
391
392 //=======================================================================
393 //function : writePoints
394 //purpose  : 
395 //=======================================================================
396
397 static bool writePoints (ofstream &                       theFile,
398                          SMESHDS_Mesh *                   theMesh,
399                          map <int,int> &                  theSmdsToGhs3dIdMap,
400                          map <int,const SMDS_MeshNode*> & theGhs3dIdToNodeMap)
401 {
402   // record structure:
403   //
404   // NB_NODES
405   // Loop from 1 to NB_NODES
406   //   X Y Z DUMMY_INT
407
408   int nbNodes = theMesh->NbNodes();
409   if ( nbNodes == 0 )
410     return false;
411
412   const char* space    = "  ";
413   const int   dummyint = 0;
414
415   int aGhs3dID = 1;
416   SMDS_NodeIteratorPtr it = theMesh->nodesIterator();
417   const SMDS_MeshNode* node;
418
419   // NB_NODES
420   theFile << space << nbNodes << endl;
421   cout << endl;
422   cout << "The initial 2D mesh contains :" << endl;
423   cout << "    " << nbNodes << " vertices" << endl;
424
425   // Loop from 1 to NB_NODES
426
427   while ( it->more() )
428   {
429     node = it->next();
430     theSmdsToGhs3dIdMap.insert( map <int,int>::value_type( node->GetID(), aGhs3dID ));
431     theGhs3dIdToNodeMap.insert (map <int,const SMDS_MeshNode*>::value_type( aGhs3dID, node ));
432     aGhs3dID++;
433
434     // X Y Z DUMMY_INT
435     theFile
436       << space << node->X()
437       << space << node->Y()
438       << space << node->Z()
439       << space << dummyint;
440
441     theFile << endl;
442   }
443
444   return true;
445 }
446
447 //=======================================================================
448 //function : writePoints
449 //purpose  : 
450 //=======================================================================
451
452 static bool writePoints (ofstream &                            theFile,
453                          SMESHDS_Mesh *                        theMesh,
454                          const vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
455 {
456   // record structure:
457   //
458   // NB_NODES
459   // Loop from 1 to NB_NODES
460   //   X Y Z DUMMY_INT
461
462   //int nbNodes = theMesh->NbNodes();
463   int nbNodes = theNodeByGhs3dId.size();
464   if ( nbNodes == 0 )
465     return false;
466
467   const char* space    = "  ";
468   const int   dummyint = 0;
469
470   const SMDS_MeshNode* node;
471
472   // NB_NODES
473   theFile << space << nbNodes << endl;
474   cout << nbNodes << " nodes" << endl;
475
476   // Loop from 1 to NB_NODES
477
478   vector<const SMDS_MeshNode*>::const_iterator nodeIt = theNodeByGhs3dId.begin();
479   vector<const SMDS_MeshNode*>::const_iterator after  = theNodeByGhs3dId.end();
480   for ( ; nodeIt != after; ++nodeIt )
481   {
482     node = *nodeIt;
483
484     // X Y Z DUMMY_INT
485     theFile
486       << space << node->X()
487       << space << node->Y()
488       << space << node->Z()
489       << space << dummyint;
490
491     theFile << endl;
492   }
493
494   return true;
495 }
496
497 //=======================================================================
498 //function : readResultFile
499 //purpose  : 
500 //=======================================================================
501
502 static bool readResultFile(const int                       fileOpen,
503                            SMESHDS_Mesh*                   theMeshDS,
504                            TopoDS_Shape                    tabShape[],
505                            double**                        tabBox,
506                            const int                       nbShape,
507                            map <int,const SMDS_MeshNode*>& theGhs3dIdToNodeMap,
508                            bool                            toMeshHoles) {
509
510   struct stat status;
511   size_t      length;
512
513   char *ptr, *mapPtr;
514   char *tetraPtr;
515   char *shapePtr;
516
517   int fileStat;
518   int nbElems, nbNodes, nbInputNodes;
519   int nodeId, triangleId;
520   int nbTriangle;
521   int ID, shapeID, ghs3dShapeID;
522   int IdShapeRef = 1;
523   int compoundID =
524     nbShape ? theMeshDS->ShapeToIndex( tabShape[0] ) : theMeshDS->ShapeToIndex( theMeshDS->ShapeToMesh() );
525
526   int *tab, *tabID, *nodeID, *nodeAssigne;
527   double *coord;
528   const SMDS_MeshNode **node;
529
530   tab    = new int[3];
531   //tabID  = new int[nbShape];
532   nodeID = new int[4];
533   coord  = new double[3];
534   node   = new const SMDS_MeshNode*[4];
535
536   TopoDS_Shape aSolid;
537   SMDS_MeshNode * aNewNode;
538   map <int,const SMDS_MeshNode*>::iterator itOnNode;
539   SMDS_MeshElement* aTet;
540 #ifdef _DEBUG_
541   set<int> shapeIDs;
542 #endif
543
544   // Read the file state
545   fileStat = fstat(fileOpen, &status);
546   length   = status.st_size;
547
548   // Mapping the result file into memory
549   ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
550   mapPtr = ptr;
551
552   ptr      = readMapIntLine(ptr, tab);
553   tetraPtr = ptr;
554
555   nbElems      = tab[0];
556   nbNodes      = tab[1];
557   nbInputNodes = tab[2];
558
559   nodeAssigne = new int[ nbNodes+1 ];
560
561   if (nbShape > 0)
562     aSolid = tabShape[0];
563
564   // Reading the nodeId
565   for (int i=0; i < 4*nbElems; i++)
566     nodeId = strtol(ptr, &ptr, 10);
567
568   // Reading the nodeCoor and update the nodeMap
569   for (int iNode=1; iNode <= nbNodes; iNode++) {
570     for (int iCoor=0; iCoor < 3; iCoor++)
571       coord[ iCoor ] = strtod(ptr, &ptr);
572     nodeAssigne[ iNode ] = 1;
573     if ( iNode > nbInputNodes ) {
574       nodeAssigne[ iNode ] = 0;
575       aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
576       theGhs3dIdToNodeMap.insert(theGhs3dIdToNodeMap.end(), make_pair( iNode, aNewNode ));
577     }
578   }
579
580   // Reading the number of triangles which corresponds to the number of shapes
581   nbTriangle = strtol(ptr, &ptr, 10);
582
583   tabID  = new int[nbTriangle];
584   for (int i=0; i<nbTriangle; i++)
585     tabID[i] = 0;
586
587   for (int i=0; i < 3*nbTriangle; i++)
588     triangleId = strtol(ptr, &ptr, 10);
589
590   shapePtr = ptr;
591
592   if ( nbTriangle <= nbShape ) // no holes
593     toMeshHoles = true; // not avoid creating tetras in holes
594
595   // Associating the tetrahedrons to the shapes
596   shapeID = compoundID;
597   for (int iElem = 0; iElem < nbElems; iElem++) {
598     for (int iNode = 0; iNode < 4; iNode++) {
599       ID = strtol(tetraPtr, &tetraPtr, 10);
600       itOnNode = theGhs3dIdToNodeMap.find(ID);
601       node[ iNode ] = itOnNode->second;
602       nodeID[ iNode ] = ID;
603     }
604     // We always run GHS3D with "to mesh holes'==TRUE but we must not create
605     // tetras within holes depending on hypo option,
606     // so we first check if aTet is inside a hole and then create it 
607     //aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
608     if ( nbTriangle > 1 ) {
609       shapeID = -1; // negative shapeID means not to create tetras if !toMeshHoles
610       ghs3dShapeID = strtol(shapePtr, &shapePtr, 10) - IdShapeRef;
611       if ( tabID[ ghs3dShapeID ] == 0 ) {
612         TopAbs_State state;
613         aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape, &state);
614         if ( toMeshHoles || state == TopAbs_IN )
615           shapeID = theMeshDS->ShapeToIndex( aSolid );
616         tabID[ ghs3dShapeID ] = shapeID;
617       }
618       else
619         shapeID = tabID[ ghs3dShapeID ];
620     }
621     else if ( nbShape > 1 ) {
622       // Case where nbTriangle == 1 while nbShape == 2 encountered
623       // with compound of 2 boxes and "To mesh holes"==False,
624       // so there are no subdomains specified for each tetrahedron.
625       // Try to guess a solid by a node already bound to shape
626       shapeID = 0;
627       for ( int i=0; i<4 && shapeID==0; i++ ) {
628         if ( nodeAssigne[ nodeID[i] ] == 1 &&
629              node[i]->GetPosition()->GetTypeOfPosition() == SMDS_TOP_3DSPACE &&
630              node[i]->GetPosition()->GetShapeId() > 1 )
631         {
632           shapeID = node[i]->GetPosition()->GetShapeId();
633         }
634       }
635       if ( shapeID==0 ) {
636         aSolid = findShape(node, aSolid, tabShape, tabBox, nbShape);
637         shapeID = theMeshDS->ShapeToIndex( aSolid );
638       }
639     }
640     // set new nodes and tetrahedron onto the shape
641     for ( int i=0; i<4; i++ ) {
642       if ( nodeAssigne[ nodeID[i] ] == 0 ) {
643         if ( shapeID > 0 )
644           theMeshDS->SetNodeInVolume( node[i], shapeID );
645         nodeAssigne[ nodeID[i] ] = shapeID;
646       }
647     }
648     if ( toMeshHoles || shapeID > 0 ) {
649       aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
650       theMeshDS->SetMeshElementOnShape( aTet, shapeID );
651     }
652 #ifdef _DEBUG_
653     shapeIDs.insert( shapeID );
654 #endif
655   }
656   // Remove nodes of tetras inside holes if !toMeshHoles
657   if ( !toMeshHoles ) {
658     itOnNode = theGhs3dIdToNodeMap.find( nbInputNodes );
659     for ( ; itOnNode != theGhs3dIdToNodeMap.end(); ++itOnNode) {
660       ID = itOnNode->first;
661       if ( nodeAssigne[ ID ] < 0 )
662         theMeshDS->RemoveFreeNode( itOnNode->second, 0 );
663     }
664   }
665
666   if ( nbElems )
667     cout << nbElems << " tetrahedrons have been associated to " << nbShape << " shapes" << endl;
668   munmap(mapPtr, length);
669   close(fileOpen);
670
671   delete [] tab;
672   delete [] tabID;
673   delete [] nodeID;
674   delete [] coord;
675   delete [] node;
676   delete [] nodeAssigne;
677
678 #ifdef _DEBUG_
679   if ( shapeIDs.size() != nbShape ) {
680     cout << "Only " << shapeIDs.size() << " solids of " << nbShape << " found" << endl;
681     for (int i=0; i<nbShape; i++) {
682       shapeID = theMeshDS->ShapeToIndex( tabShape[i] );
683       if ( shapeIDs.find( shapeID ) == shapeIDs.end() )
684         cout << "  Solid #" << shapeID << " not found" << endl;
685     }
686   }
687 #endif
688
689   return true;
690 }
691
692 //=======================================================================
693 //function : readResultFile
694 //purpose  : 
695 //=======================================================================
696
697 static bool readResultFile(const int                      fileOpen,
698                            SMESHDS_Mesh*                  theMeshDS,
699                            TopoDS_Shape                   aSolid,
700                            vector <const SMDS_MeshNode*>& theNodeByGhs3dId) {
701
702   struct stat  status;
703   size_t       length;
704
705   char *ptr, *mapPtr;
706   char *tetraPtr;
707   char *shapePtr;
708
709   int fileStat;
710   int nbElems, nbNodes, nbInputNodes;
711   int nodeId, triangleId;
712   int nbTriangle;
713   int ID, shapeID;
714
715   int *tab;
716   double *coord;
717   const SMDS_MeshNode **node;
718
719   tab   = new int[3];
720   coord = new double[3];
721   node  = new const SMDS_MeshNode*[4];
722
723   SMDS_MeshNode * aNewNode;
724   map <int,const SMDS_MeshNode*>::iterator IdNode;
725   SMDS_MeshElement* aTet;
726
727   // Read the file state
728   fileStat = fstat(fileOpen, &status);
729   length   = status.st_size;
730
731   // Mapping the result file into memory
732   ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
733   mapPtr = ptr;
734
735   ptr      = readMapIntLine(ptr, tab);
736   tetraPtr = ptr;
737
738   nbElems      = tab[0];
739   nbNodes      = tab[1];
740   nbInputNodes = tab[2];
741
742   theNodeByGhs3dId.resize( nbNodes );
743
744   // Reading the nodeId
745   for (int i=0; i < 4*nbElems; i++)
746     nodeId = strtol(ptr, &ptr, 10);
747
748   // Reading the nodeCoor and update the nodeMap
749   shapeID = theMeshDS->ShapeToIndex( aSolid );
750   for (int iNode=0; iNode < nbNodes; iNode++) {
751     for (int iCoor=0; iCoor < 3; iCoor++)
752       coord[ iCoor ] = strtod(ptr, &ptr);
753     if ((iNode+1) > nbInputNodes) {
754       aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
755       theMeshDS->SetNodeInVolume( aNewNode, shapeID );
756       theNodeByGhs3dId[ iNode ] = aNewNode;
757     }
758   }
759
760   // Reading the triangles
761   nbTriangle = strtol(ptr, &ptr, 10);
762
763   for (int i=0; i < 3*nbTriangle; i++)
764     triangleId = strtol(ptr, &ptr, 10);
765
766   shapePtr = ptr;
767
768   // Associating the tetrahedrons to the shapes
769   for (int iElem = 0; iElem < nbElems; iElem++) {
770     for (int iNode = 0; iNode < 4; iNode++) {
771       ID = strtol(tetraPtr, &tetraPtr, 10);
772       node[ iNode ] = theNodeByGhs3dId[ ID-1 ];
773     }
774     aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
775     shapeID = theMeshDS->ShapeToIndex( aSolid );
776     theMeshDS->SetMeshElementOnShape( aTet, shapeID );
777   }
778   if ( nbElems )
779     cout << nbElems << " tetrahedrons have been associated to " << nbTriangle << " shapes" << endl;
780   munmap(mapPtr, length);
781   close(fileOpen);
782
783   delete [] tab;
784   delete [] coord;
785   delete [] node;
786
787   return true;
788 }
789
790 //================================================================================
791 /*!
792  * \brief Look for a line containing a text in a file
793   * \retval bool - true if the line is found
794  */
795 //================================================================================
796
797 static bool findLineContaing(const TCollection_AsciiString& theText,
798                              const TCollection_AsciiString& theFile,
799                              TCollection_AsciiString &      theFoundLine)
800 {
801   bool found = false;
802   if ( FILE * aFile = fopen( theFile.ToCString(), "r" ))
803   {
804     char * aPtr;
805     char aBuffer[ GHS3DPlugin_BUFLENGTH ];
806     int aLineNb = 0;
807     do {
808       GHS3DPlugin_ReadLine( aPtr, aBuffer, aFile, aLineNb );
809       if ( aPtr ) {
810         theFoundLine = aPtr;
811         found = theFoundLine.Search( theText ) >= 0;
812       }
813     } while ( aPtr && !found );
814
815     fclose( aFile );
816   }
817   return found;
818 }
819
820 //================================================================================
821 /*!
822  * \brief Provide human readable text by error code reported by ghs3d
823  */
824 //================================================================================
825
826 static TCollection_AsciiString translateError(const TCollection_AsciiString& errorLine)
827 {
828   int beg = errorLine.Location("ERR ", 1, errorLine.Length());
829   if ( !beg ) return errorLine;
830
831   TCollection_AsciiString errCodeStr = errorLine.SubString( beg + 4 , errorLine.Length());
832   if ( !errCodeStr.IsIntegerValue() )
833     return errorLine;
834   Standard_Integer errCode = errCodeStr.IntegerValue();
835
836   int codeEnd = beg + 7;
837   while ( codeEnd < errorLine.Length() && isdigit( errorLine.Value(codeEnd+1) ))
838     codeEnd++;
839   TCollection_AsciiString error = errorLine.SubString( beg, codeEnd ) + ": ";
840
841   switch ( errCode ) {
842   case 0:
843     return "The surface mesh includes a face of type type other than edge, "
844       "triangle or quadrilateral. This face type is not supported.";
845   case 1:
846     return error + "Not enough memory for the face table";
847   case 2:
848     return error + "Not enough memory";
849   case 3:
850     return error + "Not enough memory";
851   case 4:
852     return error + "Face is ignored";
853   case 5:
854     return error + errorLine.SubString( beg + 5, errorLine.Length()) +
855       ": End of file. Some data are missing in the file.";
856   case 6:
857     return error + errorLine.SubString( beg + 5, errorLine.Length()) +
858       ": Read error on the file. There are wrong data in the file.";
859   case 7:
860     return error + "the metric file is inadequate (dimension other than 3).";
861   case 8:
862     return error + "the metric file is inadequate (values not per vertices).";
863   case 9:
864       return error + "the metric file contains more than one field.";
865   case 10:
866     return error + "the number of values in the \".bb\" (metric file) is incompatible with the expected"
867       "value of number of mesh vertices in the \".noboite\" file.";
868   case 12:
869     return error + "Too many sub-domains";
870   case 13:
871     return error + "the number of vertices is negative or null.";
872   case 14:
873     return error + "the number of faces is negative or null.";
874   case 15:
875     return error + "A facehas a null vertex.";
876   case 22:
877     return error + "incompatible data";
878   case 131:
879     return error + "the number of vertices is negative or null.";
880   case 132:
881     return error + "the number of vertices is negative or null (in the \".mesh\" file).";
882   case 133:
883     return error + "the number of faces is negative or null.";
884   case 1000:
885     return error + "A face appears more than once in the input surface mesh.";
886   case 1001:
887     return error + "An edge appears more than once in the input surface mesh.";
888   case 1002:
889     return error + "A face has a vertex negative or null.";
890   case 1003:
891     return error + "NOT ENOUGH MEMORY";
892   case 2000:
893     return error + "Not enough available memory.";
894   case 2002:
895     return error + "Some initial points cannot be inserted. The surface mesh is probably very bad "
896       "in terms of quality or the input list of points is wrong";
897   case 2003:
898     return error + "Some vertices are too close to one another or coincident.";
899   case 2004:
900     return error + "Some vertices are too close to one another or coincident.";
901   case 2012:
902     return error + "A vertex cannot be inserted.";
903   case 2014:
904     return error + "There are at least two points considered as coincident";
905   case 2103:
906     return error + "Some vertices are too close to one another or coincident";
907   case 3000:
908     return error + "The surface mesh regeneration step has failed";
909   case 3009:
910     return error + "Constrained edge cannot be enforced";
911   case 3019:
912     return error + "Constrained face cannot be enforced";
913   case 3029:
914     return error + "Missing faces";
915   case 3100:
916     return error + "No guess to start the definition of the connected component(s)";
917   case 3101:
918     return error + "The surface mesh includes at least one hole. The domain is not well defined";
919   case 3102:
920     return error + "Impossible to define a component";
921   case 3103:
922     return error + "The surface edge intersects another surface edge";
923   case 3104:
924     return error + "The surface edge intersects the surface face";
925   case 3105:
926     return error + "One boundary point lies within a surface face";
927   case 3106:
928     return error + "One surface edge intersects a surface face";
929   case 3107:
930     return error + "One boundary point lies within a surface edge";
931   case 3108:
932     return error + "Insufficient memory ressources detected due to a bad quality surface mesh leading "
933       "to too many swaps";
934   case 3109:
935     return error + "Edge is unique (i.e., bounds a hole in the surface)";
936   case 3122:
937     return error + "Presumably, the surface mesh is not compatible with the domain being processed";
938   case 3123:
939     return error + "Too many components, too many sub-domain";
940   case 3209:
941     return error + "The surface mesh includes at least one hole. "
942       "Therefore there is no domain properly defined";
943   case 3300:
944     return error + "Statistics.";
945   case 3400:
946     return error + "Statistics.";
947   case 3500:
948     return error + "Warning, it is dramatically tedious to enforce the boundary items";
949   case 4000:
950     return error + "Not enough memory at this time, nevertheless, the program continues. "
951       "The expected mesh will be correct but not really as large as required";
952   case 4002:
953     return error + "see above error code, resulting quality may be poor.";
954   case 4003:
955     return error + "Not enough memory at this time, nevertheless, the program continues (warning)";
956   case 8000:
957     return error + "Unknown face type.";
958   case 8005:
959   case 8006:
960     return error + errorLine.SubString( beg + 5, errorLine.Length()) +
961       ": End of file. Some data are missing in the file.";
962   case 9000:
963     return error + "A too small volume element is detected";
964   case 9001:
965     return error + "There exists at least a null or negative volume element";
966   case 9002:
967     return error + "There exist null or negative volume elements";
968   case 9003:
969     return error + "A too small volume element is detected. A face is considered being degenerated";
970   case 9100:
971     return error + "Some element is suspected to be very bad shaped or wrong";
972   case 9102:
973     return error + "A too bad quality face is detected. This face is considered degenerated";
974   case 9112:
975     return error + "A too bad quality face is detected. This face is degenerated";
976   case 9122:
977     return error + "Presumably, the surface mesh is not compatible with the domain being processed";
978   case 9999:
979     return error + "Abnormal error occured, contact hotline";
980   case 23600:
981     return error + "Not enough memory for the face table";
982   case 23601:
983     return error + "The algorithm cannot run further. "
984       "The surface mesh is probably very bad in terms of quality.";
985   case 23602:
986     return error + "Bad vertex number";
987   }
988   return errorLine;
989 }
990
991 //=============================================================================
992 /*!
993  *Here we are going to use the GHS3D mesher
994  */
995 //=============================================================================
996
997 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh&         theMesh,
998                                 const TopoDS_Shape& theShape)
999 {
1000   bool Ok(false);
1001   SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1002
1003   _nbShape = countShape( meshDS, TopAbs_SOLID ); // we count the number of shapes
1004
1005   // create bounding box for every shape inside the compound
1006
1007   int iShape = 0;
1008   TopoDS_Shape* tabShape;
1009   double**      tabBox;
1010   tabShape = new TopoDS_Shape[_nbShape];
1011   tabBox   = new double*[_nbShape];
1012   for (int i=0; i<_nbShape; i++)
1013     tabBox[i] = new double[6];
1014   Standard_Real Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
1015
1016   TopExp_Explorer expBox (meshDS->ShapeToMesh(), TopAbs_SOLID);
1017   for (; expBox.More(); expBox.Next()) {
1018     tabShape[iShape] = expBox.Current();
1019     Bnd_Box BoundingBox;
1020     BRepBndLib::Add(expBox.Current(), BoundingBox);
1021     BoundingBox.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
1022     tabBox[iShape][0] = Xmin; tabBox[iShape][1] = Xmax;
1023     tabBox[iShape][2] = Ymin; tabBox[iShape][3] = Ymax;
1024     tabBox[iShape][4] = Zmin; tabBox[iShape][5] = Zmax;
1025     iShape++;
1026   }
1027
1028   // a unique working file name
1029   // to avoid access to the same files by eg different users
1030   TCollection_AsciiString aGenericName
1031     = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1032
1033   TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1034   TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1035   aFacesFileName  = aGenericName + ".faces";  // in faces
1036   aPointsFileName = aGenericName + ".points"; // in points
1037   aResultFileName = aGenericName + ".noboite";// out points and volumes
1038   aBadResFileName = aGenericName + ".boite";  // out bad result
1039   aBbResFileName  = aGenericName + ".bb";     // out vertex stepsize
1040   aLogFileName    = aGenericName + ".log";    // log
1041
1042   // -----------------
1043   // make input files
1044   // -----------------
1045
1046   ofstream aFacesFile  ( aFacesFileName.ToCString()  , ios::out);
1047   ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
1048
1049   Ok =
1050 #ifdef WIN32
1051     aFacesFile->is_open() && aPointsFile->is_open();
1052 #else
1053     aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1054 #endif
1055   if (!Ok) {
1056     INFOS( "Can't write into " << aFacesFileName);
1057     return error(SMESH_Comment("Can't write into ") << aFacesFileName);
1058   }
1059   map <int,int> aSmdsToGhs3dIdMap;
1060   map <int,const SMDS_MeshNode*> aGhs3dIdToNodeMap;
1061
1062   Ok = writePoints( aPointsFile, meshDS, aSmdsToGhs3dIdMap, aGhs3dIdToNodeMap ) &&
1063        writeFaces ( aFacesFile,  meshDS, aSmdsToGhs3dIdMap );
1064
1065   aFacesFile.close();
1066   aPointsFile.close();
1067
1068   if ( ! Ok ) {
1069     if ( !_keepFiles ) {
1070       OSD_File( aFacesFileName ).Remove();
1071       OSD_File( aPointsFileName ).Remove();
1072     }
1073     return error(COMPERR_BAD_INPUT_MESH);
1074   }
1075   OSD_File( aResultFileName ).Remove(); // needed for boundary recovery module usage
1076
1077   // -----------------
1078   // run ghs3d mesher
1079   // -----------------
1080
1081   TCollection_AsciiString cmd( (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp ).c_str() );
1082   cmd += TCollection_AsciiString(" -f ") + aGenericName;  // file to read
1083   cmd += TCollection_AsciiString(" 1>" ) + aLogFileName;  // dump into file
1084
1085   cout << endl;
1086   cout << "Ghs3d execution..." << endl;
1087   cout << cmd << endl;
1088
1089   system( cmd.ToCString() ); // run
1090
1091   cout << endl;
1092   cout << "End of Ghs3d execution !" << endl;
1093
1094   // --------------
1095   // read a result
1096   // --------------
1097
1098   // Mapping the result file
1099
1100   int fileOpen;
1101   fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1102   if ( fileOpen < 0 ) {
1103     cout << endl;
1104     cout << "Can't open the " << aResultFileName.ToCString() << " GHS3D output file" << endl;
1105     Ok = false;
1106   }
1107   else {
1108     bool toMeshHoles =
1109       _hyp ? _hyp->GetToMeshHoles(true) : GHS3DPlugin_Hypothesis::DefaultMeshHoles();
1110     Ok = readResultFile( fileOpen, meshDS, tabShape, tabBox, _nbShape, aGhs3dIdToNodeMap,
1111                          toMeshHoles );
1112   }
1113
1114   // ---------------------
1115   // remove working files
1116   // ---------------------
1117
1118   if ( Ok )
1119   {
1120     if ( !_keepFiles )
1121       OSD_File( aLogFileName ).Remove();
1122   }
1123   else if ( OSD_File( aLogFileName ).Size() > 0 )
1124   {
1125     INFOS( "GHS3D Error, see the " << aLogFileName.ToCString() << " file" );
1126
1127     // get problem description from the log file
1128     SMESH_Comment comment;
1129     TCollection_AsciiString foundLine;
1130     if ( findLineContaing( "has expired",aLogFileName,foundLine) &&
1131          foundLine.Search("Licence") >= 0)
1132     {
1133       foundLine.LeftAdjust();
1134       comment << foundLine;
1135     }
1136     if ( findLineContaing( "%% ERROR",aLogFileName,foundLine) ||
1137          findLineContaing( " ERR ",aLogFileName,foundLine))
1138     {
1139       foundLine.LeftAdjust();
1140       comment << translateError( foundLine );
1141     }
1142     else if ( findLineContaing( "%% NO SAVING OPERATION",aLogFileName,foundLine))
1143     {
1144       comment << "Too many elements generated for a trial version.\n";
1145     }
1146     if ( comment.empty() )
1147       comment << "See " << aLogFileName << " for problem description";
1148     else
1149       comment << "\nSee " << aLogFileName << " for more information";
1150     error(COMPERR_ALGO_FAILED, comment);
1151   }
1152   else
1153   {
1154     // the log file is empty
1155     OSD_File( aLogFileName ).Remove();
1156     INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1157     error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1158   }
1159
1160   if ( !_keepFiles ) {
1161     OSD_File( aFacesFileName ).Remove();
1162     OSD_File( aPointsFileName ).Remove();
1163     OSD_File( aResultFileName ).Remove();
1164     OSD_File( aBadResFileName ).Remove();
1165     OSD_File( aBbResFileName ).Remove();
1166   }
1167   cout << "<" << aResultFileName.ToCString() << "> GHS3D output file ";
1168   if ( !Ok )
1169     cout << "not ";
1170   cout << "treated !" << endl;
1171   cout << endl;
1172
1173   _nbShape = 0;    // re-initializing _nbShape for the next Compute() method call
1174   delete [] tabShape;
1175   delete [] tabBox;
1176
1177   return Ok;
1178 }
1179
1180 //=============================================================================
1181 /*!
1182  *Here we are going to use the GHS3D mesher w/o geometry
1183  */
1184 //=============================================================================
1185 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh&         theMesh,
1186                                 SMESH_MesherHelper* aHelper)
1187 {
1188   MESSAGE("GHS3DPlugin_GHS3D::Compute()");
1189
1190   SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
1191   TopoDS_Shape theShape = aHelper->GetSubShape();
1192
1193   // a unique working file name
1194   // to avoid access to the same files by eg different users
1195   TCollection_AsciiString aGenericName
1196     = (char*) GHS3DPlugin_Hypothesis::GetFileName(_hyp).c_str();
1197
1198   TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
1199   TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
1200   aFacesFileName  = aGenericName + ".faces";  // in faces
1201   aPointsFileName = aGenericName + ".points"; // in points
1202   aResultFileName = aGenericName + ".noboite";// out points and volumes
1203   aBadResFileName = aGenericName + ".boite";  // out bad result
1204   aBbResFileName  = aGenericName + ".bb";     // out vertex stepsize
1205   aLogFileName    = aGenericName + ".log";    // log
1206
1207   // -----------------
1208   // make input files
1209   // -----------------
1210
1211   ofstream aFacesFile  ( aFacesFileName.ToCString()  , ios::out);
1212   ofstream aPointsFile  ( aPointsFileName.ToCString()  , ios::out);
1213   bool Ok =
1214 #ifdef WIN32
1215     aFacesFile->is_open() && aPointsFile->is_open();
1216 #else
1217     aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
1218 #endif
1219
1220   if (!Ok)
1221     return error( SMESH_Comment("Can't write into ") << aPointsFileName);
1222
1223   vector <const SMDS_MeshNode*> aNodeByGhs3dId;
1224
1225   Ok = (writeFaces ( aFacesFile, meshDS, aNodeByGhs3dId ) &&
1226         writePoints( aPointsFile, meshDS, aNodeByGhs3dId));
1227   
1228   aFacesFile.close();
1229   aPointsFile.close();
1230   
1231   if ( ! Ok ) {
1232     if ( !_keepFiles ) {
1233       OSD_File( aFacesFileName ).Remove();
1234       OSD_File( aPointsFileName ).Remove();
1235     }
1236     return error(COMPERR_BAD_INPUT_MESH);
1237   }
1238   OSD_File( aResultFileName ).Remove(); // needed for boundary recovery module usage
1239
1240   // -----------------
1241   // run ghs3d mesher
1242   // -----------------
1243
1244   TCollection_AsciiString cmd =
1245     (char*)GHS3DPlugin_Hypothesis::CommandToRun( _hyp, false ).c_str();
1246   cmd += TCollection_AsciiString(" -f ") + aGenericName;  // file to read
1247   cmd += TCollection_AsciiString(" 1>" ) + aLogFileName;  // dump into file
1248   
1249   system( cmd.ToCString() ); // run
1250
1251   // --------------
1252   // read a result
1253   // --------------
1254   int fileOpen;
1255   fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
1256   if ( fileOpen < 0 ) {
1257     cout << endl;
1258     cout << "Error when opening the " << aResultFileName.ToCString() << " file" << endl;
1259     cout << endl;
1260     Ok = false;
1261   }
1262   else {
1263     Ok = readResultFile( fileOpen, meshDS, theShape ,aNodeByGhs3dId );
1264   }
1265   
1266   // ---------------------
1267   // remove working files
1268   // ---------------------
1269
1270   if ( Ok )
1271   {
1272     if ( !_keepFiles )
1273       OSD_File( aLogFileName ).Remove();
1274   }
1275   else if ( OSD_File( aLogFileName ).Size() > 0 )
1276   {
1277     INFOS( "GHS3D Error, see the " << aLogFileName.ToCString() << " file" );
1278
1279     // get problem description from the log file
1280     SMESH_Comment comment;
1281     TCollection_AsciiString foundLine;
1282     if ( findLineContaing( "has expired",aLogFileName,foundLine) &&
1283          foundLine.Search("Licence") >= 0)
1284     {
1285       foundLine.LeftAdjust();
1286       comment << foundLine;
1287     }
1288     if ( findLineContaing( "%% ERROR",aLogFileName,foundLine) ||
1289          findLineContaing( " ERR ",aLogFileName,foundLine))
1290     {
1291       foundLine.LeftAdjust();
1292       comment << translateError( foundLine );
1293     }
1294     else if ( findLineContaing( "%% NO SAVING OPERATION",aLogFileName,foundLine))
1295     {
1296       comment << "Too many elements generated for a trial version.\n";
1297     }
1298     if ( comment.empty() )
1299       comment << "See " << aLogFileName << " for problem description";
1300     else
1301       comment << "\nSee " << aLogFileName << " for more information";
1302     error(COMPERR_ALGO_FAILED, comment);
1303   }
1304   else {
1305     // the log file is empty
1306     OSD_File( aLogFileName ).Remove();
1307     INFOS( "GHS3D Error, command '" << cmd.ToCString() << "' failed" );
1308     error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
1309   }
1310
1311   if ( !_keepFiles )
1312   {
1313     OSD_File( aFacesFileName ).Remove();
1314     OSD_File( aPointsFileName ).Remove();
1315     OSD_File( aResultFileName ).Remove();
1316     OSD_File( aBadResFileName ).Remove();
1317     OSD_File( aBbResFileName ).Remove();
1318   }
1319   
1320   return Ok;
1321 }