Salome HOME
NPAL16132 (GHS3D execution failed in Salome version 4.0.0)
[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 "SMESH_Gen.hxx"
31 #include "SMESH_Mesh.hxx"
32 #include "SMESH_Comment.hxx"
33
34 #include "SMDS_MeshElement.hxx"
35 #include "SMDS_MeshNode.hxx"
36
37 #include <TopExp_Explorer.hxx>
38 #include <OSD_File.hxx>
39
40 #include "utilities.h"
41
42 #ifndef WIN32
43 #include <sys/sysinfo.h>
44 #endif
45
46 //#include <Standard_Stream.hxx>
47
48 #include <BRepGProp.hxx>
49 #include <BRepBndLib.hxx>
50 #include <BRepClass_FaceClassifier.hxx>
51 #include <BRepClass3d_SolidClassifier.hxx>
52 #include <TopAbs.hxx>
53 #include <Bnd_Box.hxx>
54 #include <GProp_GProps.hxx>
55 #include <Precision.hxx>
56
57 #define castToNode(n) static_cast<const SMDS_MeshNode *>( n );
58
59 #ifdef _DEBUG_
60 #define DUMP(txt) \
61 //  cout << txt
62 #else
63 #define DUMP(txt)
64 #endif
65
66 extern "C"
67 {
68 #ifndef WNT
69 #include <unistd.h>
70 #include <sys/mman.h>
71 #endif
72 #include <sys/stat.h>
73 #include <fcntl.h>
74 }
75
76 //=============================================================================
77 /*!
78  *  
79  */
80 //=============================================================================
81
82 GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D(int hypId, int studyId, SMESH_Gen* gen)
83   : SMESH_3D_Algo(hypId, studyId, gen)
84 {
85   MESSAGE("GHS3DPlugin_GHS3D::GHS3DPlugin_GHS3D");
86   _name = "GHS3D_3D";
87   _shapeType = (1 << TopAbs_SHELL) | (1 << TopAbs_SOLID);// 1 bit /shape type
88   _onlyUnaryInput = false; // Compute() will be called on a compound of solids
89 //   _iShape=0;
90 //   _nbShape=0;
91 }
92
93 //=============================================================================
94 /*!
95  *  
96  */
97 //=============================================================================
98
99 GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D()
100 {
101   MESSAGE("GHS3DPlugin_GHS3D::~GHS3DPlugin_GHS3D");
102 }
103
104 //=============================================================================
105 /*!
106  *  
107  */
108 //=============================================================================
109
110 bool GHS3DPlugin_GHS3D::CheckHypothesis ( SMESH_Mesh&                          aMesh,
111                                           const TopoDS_Shape&                  aShape,
112                                           SMESH_Hypothesis::Hypothesis_Status& aStatus )
113 {
114 //  MESSAGE("GHS3DPlugin_GHS3D::CheckHypothesis");
115   aStatus = SMESH_Hypothesis::HYP_OK;
116   return true;
117 }
118
119 //================================================================================
120 /*!
121  * \brief Write faces bounding theShape to file
122  */
123 //================================================================================
124
125 static bool writeFaces (ofstream &                      theFile,
126                         SMESHDS_Mesh *                  theMesh,
127                         const TopoDS_Shape&             theShape,
128                         vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
129 {
130   // record structure:
131   //
132   // NB_ELEMS DUMMY_INT
133   // Loop from 1 to NB_ELEMS
134   // NB_NODES NODE_NB_1 NODE_NB_2 ... (NB_NODES + 1) times: DUMMY_INT
135
136   // get all faces bound to theShape
137
138   // Solids in the ShapeToMesh() can be meshed by different meshers,
139   // so we take faces only from the given shape
140   //TopoDS_Shape theShape = theMesh->ShapeToMesh();
141   int nbFaces = 0;
142   list< const SMDS_MeshElement* > faces;
143   // Use TopTools_IndexedMapOfShape in order not to take twice mesh faces from
144   // a geom face shared by two solids
145   TopTools_IndexedMapOfShape faceMap;
146   TopExp::MapShapes( theShape, TopAbs_FACE, faceMap );
147   SMESHDS_SubMesh* sm;
148   SMDS_ElemIteratorPtr eIt;
149
150   const char* space    = "  ";
151   const int   dummyint = 0;
152
153   list< const SMDS_MeshElement* >::iterator f;
154   map< const SMDS_MeshNode*,int >::iterator it;
155   SMDS_ElemIteratorPtr nodeIt;
156   const SMDS_MeshElement* elem;
157   int nbNodes;
158   //int aSmdsID;
159
160   for ( int i = 0; i < faceMap.Extent(); ++i ) {
161     sm = theMesh->MeshElements( faceMap( i+1 ) );
162     if ( sm ) {
163       eIt = sm->GetElements();
164       while ( eIt->more() ) {
165         faces.push_back( eIt->next() );
166         nbFaces++;
167       }
168     }
169   }
170
171   if ( nbFaces == 0 )
172     return false;
173
174   cout << "The initial 2D mesh contains " << nbFaces << " faces and ";
175
176   // NB_ELEMS DUMMY_INT
177   theFile << space << nbFaces << space << dummyint << endl;
178
179   // Loop from 1 to NB_ELEMS
180
181   map<const SMDS_MeshNode*,int> aNodeToGhs3dIdMap;
182   f = faces.begin();
183   for ( ; f != faces.end(); ++f )
184   {
185     // NB_NODES PER FACE
186     elem = *f;
187     nbNodes = elem->NbNodes();
188     theFile << space << nbNodes;
189
190     // NODE_NB_1 NODE_NB_2 ...
191     nodeIt = elem->nodesIterator();
192     while ( nodeIt->more() )
193     {
194       // find GHS3D ID
195       const SMDS_MeshNode* node = castToNode( nodeIt->next() );
196       int newId = aNodeToGhs3dIdMap.size() + 1; // ghs3d ids count from 1
197       it = aNodeToGhs3dIdMap.insert( make_pair( node, newId )).first;
198       theFile << space << it->second;
199     }
200
201     // (NB_NODES + 1) times: DUMMY_INT
202     for ( int i=0; i<=nbNodes; i++)
203       theFile << space << dummyint;
204
205     theFile << endl;
206   }
207
208   // put nodes to theNodeByGhs3dId vector
209   theNodeByGhs3dId.resize( aNodeToGhs3dIdMap.size() );
210   map<const SMDS_MeshNode*,int>::const_iterator n2id = aNodeToGhs3dIdMap.begin();
211   for ( ; n2id != aNodeToGhs3dIdMap.end(); ++ n2id)
212   {
213     theNodeByGhs3dId[ n2id->second - 1 ] = n2id->first; // ghs3d ids count from 1
214   }
215
216   return true;
217 }
218
219 //=======================================================================
220 //function : writePoints
221 //purpose  : 
222 //=======================================================================
223
224 static bool writePoints (ofstream &                            theFile,
225                          SMESHDS_Mesh *                        theMesh,
226                          const vector <const SMDS_MeshNode*> & theNodeByGhs3dId)
227 {
228   // record structure:
229   //
230   // NB_NODES
231   // Loop from 1 to NB_NODES
232   //   X Y Z DUMMY_INT
233
234   //int nbNodes = theMesh->NbNodes();
235   int nbNodes = theNodeByGhs3dId.size();
236   if ( nbNodes == 0 )
237     return false;
238
239   const char* space    = "  ";
240   const int   dummyint = 0;
241
242   const SMDS_MeshNode* node;
243
244   // NB_NODES
245   theFile << space << nbNodes << endl;
246   cout << nbNodes << " nodes" << endl;
247
248   // Loop from 1 to NB_NODES
249
250   vector<const SMDS_MeshNode*>::const_iterator nodeIt = theNodeByGhs3dId.begin();
251   vector<const SMDS_MeshNode*>::const_iterator after  = theNodeByGhs3dId.end();
252   for ( ; nodeIt != after; ++nodeIt )
253   {
254     node = *nodeIt;
255
256     // X Y Z DUMMY_INT
257     theFile
258       << space << node->X()
259       << space << node->Y()
260       << space << node->Z()
261       << space << dummyint;
262
263     theFile << endl;
264   }
265
266   return true;
267 }
268
269 //=======================================================================
270 //function : findSolid
271 //purpose  : 
272 //=======================================================================
273
274 static TopoDS_Shape findSolid(const SMDS_MeshNode *aNode[],
275                               TopoDS_Shape        aSolid,
276                               const TopoDS_Shape  shape[],
277                               const double        box[][6],
278                               const int           nShape) {
279
280   Standard_Real PX, PY, PZ;
281   int iShape;
282
283   PX = ( aNode[0]->X() + aNode[1]->X() + aNode[2]->X() + aNode[3]->X() ) / 4.0;
284   PY = ( aNode[0]->Y() + aNode[1]->Y() + aNode[2]->Y() + aNode[3]->Y() ) / 4.0;
285   PZ = ( aNode[0]->Z() + aNode[1]->Z() + aNode[2]->Z() + aNode[3]->Z() ) / 4.0;
286   gp_Pnt aPnt(PX, PY, PZ);
287
288   BRepClass3d_SolidClassifier SC (aSolid, aPnt, Precision::Confusion());
289   if ( not(SC.State() == TopAbs_IN) ) {
290     for (iShape = 0; iShape < nShape; iShape++) {
291       aSolid = shape[iShape];
292       if ( not( PX < box[iShape][0] || box[iShape][1] < PX ||
293                 PY < box[iShape][2] || box[iShape][3] < PY ||
294                 PZ < box[iShape][4] || box[iShape][5] < PZ) ) {
295         BRepClass3d_SolidClassifier SC (aSolid, aPnt, Precision::Confusion());
296         if (SC.State() == TopAbs_IN)
297           break;
298       }
299     }
300   }
301   return aSolid;
302 }
303
304 //=======================================================================
305 //function : readMapIntLine
306 //purpose  : 
307 //=======================================================================
308
309 static char* readMapIntLine(char* ptr, int tab[]) {
310   long int intVal;
311   cout << endl;
312
313   for ( int i=0; i<17; i++ ) {
314     intVal = strtol(ptr, &ptr, 10);
315     if ( i < 3 )
316       tab[i] = intVal;
317   }
318   return ptr;
319 }
320
321 //=======================================================================
322 //function : readLine
323 //purpose  : 
324 //=======================================================================
325
326 #define GHS3DPlugin_BUFLENGTH 256
327 #define GHS3DPlugin_ReadLine(aPtr,aBuf,aFile,aLineNb) \
328 {  aPtr = fgets( aBuf, GHS3DPlugin_BUFLENGTH - 2, aFile ); aLineNb++; DUMP(endl); }
329
330 #include <list>
331 //=======================================================================
332 //function : readResultFile
333 //purpose  : 
334 //=======================================================================
335
336 static bool readResultFile(const int                      fileOpen,
337                            SMESHDS_Mesh*                  theMeshDS,
338                            TopoDS_Shape                   tabShape[],
339                            double                         tabBox[][6],
340                            const int                      nShape,
341                            vector <const SMDS_MeshNode*>& theNodeByGhs3dId) {
342
343   struct stat  status;
344   size_t       length;
345
346   char *ptr, *mapPtr;
347   char *tetraPtr;
348   char *shapePtr;
349
350   int fileStat;
351   int nbElems, nbNodes, nbInputNodes;
352   int nodeId, triangleId;
353   int tab[3]/*, tabID[nShape]*/;
354   int nbTriangle;
355   int ID, shapeID, ghs3dShapeID;
356
357   double coord [3];
358   vector< int > tabID (nShape, 0);
359
360   TopoDS_Shape aSolid;
361   SMDS_MeshNode * aNewNode;
362   const SMDS_MeshNode * node[4];
363   map <int,const SMDS_MeshNode*>::iterator IdNode;
364   SMDS_MeshElement* aTet;
365
366 //   for (int i=0; i<nShape; i++)
367 //     tabID[i] = 0;
368
369   // Read the file state
370   fileStat = fstat(fileOpen, &status);
371   length   = status.st_size;
372
373   // Mapping the result file into memory
374   ptr = (char *) mmap(0,length,PROT_READ,MAP_PRIVATE,fileOpen,0);
375   mapPtr = ptr;
376
377   ptr      = readMapIntLine(ptr, tab);
378   tetraPtr = ptr;
379
380   nbElems      = tab[0];
381   nbNodes      = tab[1];
382   nbInputNodes = tab[2];
383
384   theNodeByGhs3dId.resize( nbNodes );
385
386   // Reading the nodeId
387   for (int i=0; i < 4*nbElems; i++)
388     nodeId = strtol(ptr, &ptr, 10);
389
390   // Reading the nodeCoor and update the nodeMap
391   for (int iNode=0; iNode < nbNodes; iNode++) {
392     for (int iCoor=0; iCoor < 3; iCoor++)
393       coord[ iCoor ] = strtod(ptr, &ptr);
394     if ((iNode+1) > nbInputNodes) {
395       aNewNode = theMeshDS->AddNode( coord[0],coord[1],coord[2] );
396       theNodeByGhs3dId[ iNode ] = aNewNode;
397     }
398   }
399
400   // Reading the triangles
401   nbTriangle = strtol(ptr, &ptr, 10);
402
403   for (int i=0; i < 3*nbTriangle; i++)
404     triangleId = strtol(ptr, &ptr, 10);
405
406   shapePtr = ptr;
407
408   // Associating the tetrahedrons to the shapes
409   for (int iElem = 0; iElem < nbElems; iElem++) {
410     for (int iNode = 0; iNode < 4; iNode++) {
411       ID = strtol(tetraPtr, &tetraPtr, 10);
412       node[ iNode ] = theNodeByGhs3dId[ ID-1 ];
413     }
414     aTet = theMeshDS->AddVolume( node[1], node[0], node[2], node[3] );
415     ghs3dShapeID = strtol(shapePtr, &shapePtr, 10);
416     if ( !ghs3dShapeID ) ghs3dShapeID = 1;
417     if ( tabID[ ghs3dShapeID - 1 ] == 0 ) {
418       if (iElem == 0)
419         aSolid = tabShape[0];
420       aSolid = findSolid(node, aSolid, tabShape, tabBox, nShape /*nbTriangle*/);
421       shapeID = theMeshDS->ShapeToIndex( aSolid );
422       tabID[ ghs3dShapeID - 1] = shapeID;
423     }
424     else {
425       shapeID = tabID[ ghs3dShapeID - 1];
426     }
427     theMeshDS->SetMeshElementOnShape( aTet, shapeID );
428     // set new nodes on to the shape
429     SMDS_ElemIteratorPtr nodeIt = aTet->nodesIterator();
430     while ( nodeIt->more() ) {
431       const SMDS_MeshNode * n = castToNode( nodeIt->next() );
432       if ( !n->GetPosition()->GetShapeId() )
433         theMeshDS->SetNodeInVolume( n, shapeID );
434     }
435   }
436   if ( nbElems )
437     cout << nbElems << " tetrahedrons have been associated to " << nbTriangle << " shapes" << endl;
438   munmap(mapPtr, length);
439   close(fileOpen);
440   return true;
441 }
442
443 //=======================================================================
444 //function : getTmpDir
445 //purpose  : 
446 //=======================================================================
447
448 static TCollection_AsciiString getTmpDir()
449 {
450   TCollection_AsciiString aTmpDir;
451
452   char *Tmp_dir = getenv("SALOME_TMP_DIR");
453   if(Tmp_dir != NULL) {
454     aTmpDir = Tmp_dir;
455 #ifdef WIN32
456     if(aTmpDir.Value(aTmpDir.Length()) != '\\') aTmpDir+='\\';
457 #else
458     if(aTmpDir.Value(aTmpDir.Length()) != '/') aTmpDir+='/';
459 #endif      
460   }
461   else {
462 #ifdef WIN32
463     aTmpDir = TCollection_AsciiString("C:\\");
464 #else
465     aTmpDir = TCollection_AsciiString("/tmp/");
466 #endif
467   }
468   return aTmpDir;
469 }
470
471 //================================================================================
472 /*!
473  * \brief Look for a line containing a text in a file
474   * \retval bool - true if the line is found
475  */
476 //================================================================================
477
478 static bool findLineContaing(const TCollection_AsciiString& theText,
479                              const TCollection_AsciiString& theFile,
480                              TCollection_AsciiString &      theFoundLine)
481 {
482   bool found = false;
483   if ( FILE * aFile = fopen( theFile.ToCString(), "r" ))
484   {
485     char * aPtr;
486     char aBuffer[ GHS3DPlugin_BUFLENGTH ];
487     int aLineNb = 0;
488     do {
489       GHS3DPlugin_ReadLine( aPtr, aBuffer, aFile, aLineNb );
490       if ( aPtr ) {
491         theFoundLine = aPtr;
492         found = theFoundLine.Search( theText ) >= 0;
493       }
494     } while ( aPtr && !found );
495
496     fclose( aFile );
497   }
498   return found;
499 }
500
501 //=============================================================================
502 /*!
503  *Here we are going to use the GHS3D mesher
504  */
505 //=============================================================================
506
507 bool GHS3DPlugin_GHS3D::Compute(SMESH_Mesh&         theMesh,
508                                 const TopoDS_Shape& theShape)
509 {
510   // theShape is a compound of solids as _onlyUnaryInput = false
511   bool Ok(false);
512   SMESHDS_Mesh* meshDS = theMesh.GetMeshDS();
513
514   int _nbShape = 0;
515   /*if (_iShape == 0 && _nbShape == 0)*/ {
516     cout << endl;
517     cout << "Ghs3d execution..." << endl;
518     cout << endl;
519     
520     //TopExp_Explorer exp (meshDS->ShapeToMesh(), TopAbs_SOLID);
521     TopExp_Explorer exp (theShape, TopAbs_SOLID);
522     for (; exp.More(); exp.Next())
523       _nbShape++;
524   }
525   
526   //_iShape++;
527
528   /*if ( _iShape == _nbShape )*/ {
529
530     // create bounding box for every shape
531
532     int iShape = 0;
533     TopoDS_Shape tabShape[_nbShape];
534     double tabBox[_nbShape][6];
535     Standard_Real Xmin, Ymin, Zmin, Xmax, Ymax, Zmax;
536
537     TopExp_Explorer expBox (theShape, TopAbs_SOLID);
538     for (; expBox.More(); expBox.Next()) {
539       tabShape[iShape] = expBox.Current();
540       Bnd_Box BoundingBox;
541       BRepBndLib::Add(expBox.Current(), BoundingBox);
542       BoundingBox.Get(Xmin, Ymin, Zmin, Xmax, Ymax, Zmax);
543       tabBox[iShape][0] = Xmin; tabBox[iShape][1] = Xmax;
544       tabBox[iShape][2] = Ymin; tabBox[iShape][3] = Ymax;
545       tabBox[iShape][4] = Zmin; tabBox[iShape][5] = Zmax;
546       iShape++;
547     }
548
549     // make a unique working file name
550     // to avoid access to the same files by eg different users
551   
552     TCollection_AsciiString aGenericName, aTmpDir = getTmpDir();
553     aGenericName = aTmpDir + "GHS3D_";
554 #ifdef WIN32
555     aGenericName += GetCurrentProcessId();
556 #else
557     aGenericName += getpid();
558 #endif
559     aGenericName += "_";
560     aGenericName += meshDS->ShapeToIndex( theShape );
561
562     TCollection_AsciiString aFacesFileName, aPointsFileName, aResultFileName;
563     TCollection_AsciiString aBadResFileName, aBbResFileName, aLogFileName;
564     aFacesFileName  = aGenericName + ".faces";  // in faces
565     aPointsFileName = aGenericName + ".points"; // in points
566     aResultFileName = aGenericName + ".noboite";// out points and volumes
567     aBadResFileName = aGenericName + ".boite";  // out bad result
568     aBbResFileName  = aGenericName + ".bb";     // out vertex stepsize
569     aLogFileName    = aGenericName + ".log";    // log
570
571     // -----------------
572     // make input files
573     // -----------------
574
575     ofstream aFacesFile  ( aFacesFileName.ToCString()  , ios::out);
576     ofstream aPointsFile ( aPointsFileName.ToCString() , ios::out);
577     // bool Ok =
578     Ok =
579 #ifdef WIN32
580       aFacesFile->is_open() && aPointsFile->is_open();
581 #else
582       aFacesFile.rdbuf()->is_open() && aPointsFile.rdbuf()->is_open();
583 #endif
584     if (!Ok)
585         return error(dfltErr(), SMESH_Comment("Can't write into ") << aTmpDir);
586
587     vector <const SMDS_MeshNode*> aNodeByGhs3dId;
588
589     Ok = ( writeFaces ( aFacesFile,  meshDS, theShape, aNodeByGhs3dId ) &&
590            writePoints( aPointsFile, meshDS, aNodeByGhs3dId ));
591
592     aFacesFile.close();
593     aPointsFile.close();
594
595     if ( ! Ok ) {
596       if ( !getenv("GHS3D_KEEP_FILES") ) {
597         OSD_File( aFacesFileName ).Remove();
598         OSD_File( aPointsFileName ).Remove();
599       }
600       return error(COMPERR_BAD_INPUT_MESH);
601     }
602
603     // -----------------
604     // run ghs3d mesher              WIN32???
605     // -----------------
606
607     // ghs3d need to know amount of memory it may use (MB).
608     // Default memory is defined at ghs3d installation but it may be not enough,
609     // so allow to use about all available memory
610
611     TCollection_AsciiString memory;
612 #ifdef WIN32
613     // ????
614 #else
615     struct sysinfo si;
616     int err = sysinfo( &si );
617     if ( !err ) {
618       int freeMem = si.totalram * si.mem_unit / 1024 / 1024;
619       memory = "-m ";
620       memory += int( 0.7 * freeMem );
621     }
622 #endif
623
624     OSD_File( aResultFileName ).Remove(); // old file prevents writing a new one
625
626     TCollection_AsciiString cmd( "ghs3d " ); // command to run
627     cmd +=
628       memory +                // memory
629       " -c 0"                 // 0 - mesh all components, 1 - keep only the main component
630       " -f " + aGenericName + // file to read
631       " 1>" + aLogFileName;   // dump into file
632
633     MESSAGE("GHS3DPlugin_GHS3D::Compute() " << cmd );
634     system( cmd.ToCString() ); // run
635
636     cout << endl;
637     cout << "End of Ghs3d execution !" << endl;
638
639     // --------------
640     // read a result
641     // --------------
642
643     // Mapping the result file
644
645     int fileOpen;
646     fileOpen = open( aResultFileName.ToCString(), O_RDONLY);
647     if ( fileOpen < 0 ) {
648       cout << endl;
649       cout << "Error when opening the " << aResultFileName.ToCString() << " file" << endl;
650       cout << endl;
651       Ok = false;
652     }
653     else {
654       Ok = readResultFile( fileOpen, meshDS, tabShape, tabBox, _nbShape, aNodeByGhs3dId );
655     }
656
657   // ---------------------
658   // remove working files
659   // ---------------------
660
661     if ( Ok )
662     {
663       OSD_File( aLogFileName ).Remove();
664     }
665     else if ( OSD_File( aLogFileName ).Size() > 0 )
666     {
667       // get problem description from the log file
668       SMESH_Comment comment;
669       TCollection_AsciiString foundLine;
670       if ( findLineContaing( "has expired",aLogFileName,foundLine) &&
671            foundLine.Search("Licence") >= 0)
672       {
673         foundLine.LeftAdjust();
674         comment << foundLine;
675       }
676       if ( findLineContaing( "%% ERROR",aLogFileName,foundLine))
677       {
678         foundLine.LeftAdjust();
679         comment << foundLine;
680       }
681       if ( findLineContaing( "%% NO SAVING OPERATION",aLogFileName,foundLine))
682       {
683         comment << "Too many elements generated for a trial version.\n";
684       }
685       if ( comment.empty() )
686         comment << "See " << aLogFileName << " for problem description";
687       else
688         comment << "See " << aLogFileName << " for more information";
689       error(COMPERR_ALGO_FAILED, comment);
690     }
691     else
692     {
693       // the log file is empty
694       OSD_File( aLogFileName ).Remove();
695       error(COMPERR_ALGO_FAILED, "ghs3d: command not found" );
696     }
697
698     if ( !getenv("GHS3D_KEEP_FILES") )
699     {
700       OSD_File( aFacesFileName ).Remove();
701       OSD_File( aPointsFileName ).Remove();
702       OSD_File( aResultFileName ).Remove();
703       OSD_File( aBadResFileName ).Remove();
704       OSD_File( aBbResFileName ).Remove();
705     }
706     /*if ( _iShape == _nbShape )*/ {
707       cout << aResultFileName.ToCString() << " Output file ";
708       if ( !Ok )
709         cout << "not ";
710       cout << "treated !" << endl;
711       cout << endl;
712     }
713   }
714
715   return Ok;
716 }